xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 6080607f3a2f62b5cd8d8dcda817fda8d27c4ec4)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
305709791SSatish Balay #include <../src/mat/impls/dense/seq/dense.h>
408122e43SStefano Zampini #include <petscblaslapack.h>
534a97f8cSStefano Zampini 
63202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
75ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
8df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
9df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
10d62866d3SStefano Zampini 
11ca92afb2SStefano Zampini /* if v2 is not present, correction is done in-place */
125cbda25cSStefano Zampini PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13ca92afb2SStefano Zampini {
14ca92afb2SStefano Zampini   PetscScalar    *array;
15ca92afb2SStefano Zampini   PetscScalar    *array2;
16ca92afb2SStefano Zampini   PetscErrorCode ierr;
17ca92afb2SStefano Zampini 
18ca92afb2SStefano Zampini   PetscFunctionBegin;
19ca92afb2SStefano Zampini   if (!ctx->benign_n) PetscFunctionReturn(0);
205cbda25cSStefano Zampini   if (sol && full) {
215cbda25cSStefano Zampini     PetscInt n_I,size_schur;
225cbda25cSStefano Zampini 
235cbda25cSStefano Zampini     /* get sizes */
24282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
255cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
265cbda25cSStefano Zampini     n_I = n_I - size_schur;
275cbda25cSStefano Zampini     /* get schur sol from array */
285cbda25cSStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
295cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
305cbda25cSStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
315cbda25cSStefano Zampini     /* apply interior sol correction */
32282d6408SStefano Zampini     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
335cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
345cbda25cSStefano Zampini     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
355cbda25cSStefano Zampini   }
36ca92afb2SStefano Zampini   if (v2) {
37ca92afb2SStefano Zampini     PetscInt nl;
38ca92afb2SStefano Zampini 
39ca92afb2SStefano Zampini     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
40ca92afb2SStefano Zampini     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
41ca92afb2SStefano Zampini     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
42ca92afb2SStefano Zampini     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
43ca92afb2SStefano Zampini   } else {
44ca92afb2SStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
45ca92afb2SStefano Zampini     array2 = array;
46ca92afb2SStefano Zampini   }
47ca92afb2SStefano Zampini   if (!sol) { /* change rhs */
48ca92afb2SStefano Zampini     PetscInt n;
49ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
50ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
51ca92afb2SStefano Zampini       const PetscInt *cols;
52ca92afb2SStefano Zampini       PetscInt       nz,i;
53ca92afb2SStefano Zampini 
54ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
55ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
56ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
5722db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
5822db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
5922db5ddcSStefano Zampini #else
60ca92afb2SStefano Zampini       sum = -sum/nz;
6122db5ddcSStefano Zampini #endif
62ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
63ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
64ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
65ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
66ca92afb2SStefano Zampini     }
67ca92afb2SStefano Zampini   } else {
68ca92afb2SStefano Zampini     PetscInt n;
69ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
70ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
71ca92afb2SStefano Zampini       const PetscInt *cols;
72ca92afb2SStefano Zampini       PetscInt       nz,i;
73ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
74ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
75ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
7622db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
7722db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
7822db5ddcSStefano Zampini #else
79ca92afb2SStefano Zampini       sum = -sum/nz;
8022db5ddcSStefano Zampini #endif
81ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
82ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
83ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
84ca92afb2SStefano Zampini     }
85ca92afb2SStefano Zampini   }
86ca92afb2SStefano Zampini   if (v2) {
87ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8872b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
89ca92afb2SStefano Zampini   } else {
90ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
91ca92afb2SStefano Zampini   }
925cbda25cSStefano Zampini   if (!sol && full) {
935cbda25cSStefano Zampini     Vec      usedv;
945cbda25cSStefano Zampini     PetscInt n_I,size_schur;
955cbda25cSStefano Zampini 
965cbda25cSStefano Zampini     /* get sizes */
97282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
985cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
995cbda25cSStefano Zampini     n_I = n_I - size_schur;
1005cbda25cSStefano Zampini     /* compute schur rhs correction */
1015cbda25cSStefano Zampini     if (v2) {
1025cbda25cSStefano Zampini       usedv = v2;
1035cbda25cSStefano Zampini     } else {
1045cbda25cSStefano Zampini       usedv = v;
1055cbda25cSStefano Zampini     }
1065cbda25cSStefano Zampini     /* apply schur rhs correction */
1075cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1085cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1095cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1105cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
111282d6408SStefano Zampini     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1125cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1135cbda25cSStefano Zampini   }
114ca92afb2SStefano Zampini   PetscFunctionReturn(0);
115ca92afb2SStefano Zampini }
116ca92afb2SStefano Zampini 
117df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118d62866d3SStefano Zampini {
119df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
120683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
121d62866d3SStefano Zampini   PetscErrorCode     ierr;
122d62866d3SStefano Zampini 
123d62866d3SStefano Zampini   PetscFunctionBegin;
124d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
125683d3df6SStefano Zampini   if (full) {
126d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
127d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
128d62866d3SStefano Zampini #endif
1295cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1305cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1315cbda25cSStefano Zampini #endif
132683d3df6SStefano Zampini     copy = ctx->has_vertices;
133d4933d67SStefano Zampini   } else { /* interior solver */
1346dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
135683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1366dba178dSStefano Zampini #endif
137d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
138d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
139d4933d67SStefano Zampini #endif
140683d3df6SStefano Zampini     copy = PETSC_TRUE;
141683d3df6SStefano Zampini   }
142683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
143683d3df6SStefano Zampini   if (copy) {
144ca92afb2SStefano Zampini     PetscInt    n;
145df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
146ca92afb2SStefano Zampini 
147ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
148683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
149df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
150df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
151df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
152683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
153683d3df6SStefano Zampini 
1545cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
155683d3df6SStefano Zampini     if (transpose) {
156683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
157683d3df6SStefano Zampini     } else {
158683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
159683d3df6SStefano Zampini     }
1605cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
161683d3df6SStefano Zampini 
162683d3df6SStefano Zampini     /* get back data to caller worskpace */
163df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
164683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
165df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
166683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
167df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
168683d3df6SStefano Zampini   } else {
169ca92afb2SStefano Zampini     if (ctx->benign_n) {
1705cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
171ca92afb2SStefano Zampini       if (transpose) {
172ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
173ca92afb2SStefano Zampini       } else {
174ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
175ca92afb2SStefano Zampini       }
1765cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
177ca92afb2SStefano Zampini     } else {
178e28d306cSStefano Zampini       if (transpose) {
179e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
180e28d306cSStefano Zampini       } else {
1816816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
182e28d306cSStefano Zampini       }
183683d3df6SStefano Zampini     }
184ca92afb2SStefano Zampini   }
1855cbda25cSStefano Zampini   /* restore defaults */
1865cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1875cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1885cbda25cSStefano Zampini #endif
189d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
190d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
191d4933d67SStefano Zampini #endif
192d62866d3SStefano Zampini   PetscFunctionReturn(0);
193d62866d3SStefano Zampini }
194d62866d3SStefano Zampini 
195df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196e28d306cSStefano Zampini {
197e28d306cSStefano Zampini   PetscErrorCode   ierr;
198e28d306cSStefano Zampini 
199e28d306cSStefano Zampini   PetscFunctionBegin;
200df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
201e28d306cSStefano Zampini   PetscFunctionReturn(0);
202e28d306cSStefano Zampini }
203e28d306cSStefano Zampini 
204df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205e28d306cSStefano Zampini {
206e28d306cSStefano Zampini   PetscErrorCode   ierr;
207e28d306cSStefano Zampini 
208e28d306cSStefano Zampini   PetscFunctionBegin;
209df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
210683d3df6SStefano Zampini   PetscFunctionReturn(0);
211683d3df6SStefano Zampini }
212683d3df6SStefano Zampini 
213df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214683d3df6SStefano Zampini {
215683d3df6SStefano Zampini   PetscErrorCode   ierr;
216683d3df6SStefano Zampini 
217683d3df6SStefano Zampini   PetscFunctionBegin;
218df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
219683d3df6SStefano Zampini   PetscFunctionReturn(0);
220683d3df6SStefano Zampini }
221683d3df6SStefano Zampini 
222df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223683d3df6SStefano Zampini {
224683d3df6SStefano Zampini   PetscErrorCode   ierr;
225683d3df6SStefano Zampini 
226683d3df6SStefano Zampini   PetscFunctionBegin;
227df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
228e28d306cSStefano Zampini   PetscFunctionReturn(0);
229e28d306cSStefano Zampini }
230e28d306cSStefano Zampini 
23115579a77SStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer)
23215579a77SStefano Zampini {
23315579a77SStefano Zampini   PCBDDCReuseSolvers ctx;
23415579a77SStefano Zampini   PetscBool          iascii;
23515579a77SStefano Zampini   PetscErrorCode     ierr;
23615579a77SStefano Zampini 
23715579a77SStefano Zampini   PetscFunctionBegin;
23815579a77SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
23915579a77SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
24015579a77SStefano Zampini   if (iascii) {
24115579a77SStefano Zampini     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
24215579a77SStefano Zampini   }
24315579a77SStefano Zampini   ierr = MatView(ctx->F,viewer);CHKERRQ(ierr);
24415579a77SStefano Zampini   if (iascii) {
24515579a77SStefano Zampini     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
24615579a77SStefano Zampini   }
24715579a77SStefano Zampini   PetscFunctionReturn(0);
24815579a77SStefano Zampini }
24915579a77SStefano Zampini 
250df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
251d62866d3SStefano Zampini {
252ca92afb2SStefano Zampini   PetscInt       i;
253d62866d3SStefano Zampini   PetscErrorCode ierr;
254d62866d3SStefano Zampini 
255d62866d3SStefano Zampini   PetscFunctionBegin;
256d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
257d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
258d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
259d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
260d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
26153892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
26253892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
263d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
26453892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
26553892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
266ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
267ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
268ca92afb2SStefano Zampini   }
269ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
270ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2715cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2725cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2735cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2745cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
275d62866d3SStefano Zampini   PetscFunctionReturn(0);
276d62866d3SStefano Zampini }
277d5574798SStefano Zampini 
2785ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2793202ece2SStefano Zampini {
2803202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2813202ece2SStefano Zampini   KSP            ksp;
2823202ece2SStefano Zampini   PC             pc;
2833202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2843202ece2SStefano Zampini   PetscReal      fill = 2.0;
285f11841e3SStefano Zampini   PetscInt       n_I;
2863202ece2SStefano Zampini   PetscMPIInt    size;
2873202ece2SStefano Zampini   PetscErrorCode ierr;
2883202ece2SStefano Zampini 
2893202ece2SStefano Zampini   PetscFunctionBegin;
2903202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2916c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
292f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
293f11841e3SStefano Zampini     PetscBool Sdense;
294f11841e3SStefano Zampini 
295f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2966c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
297f11841e3SStefano Zampini   }
2983202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2993202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
3003202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
3013202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
3023202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
3033202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
3043202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
3053202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
306f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
307f11841e3SStefano Zampini   if (n_I) {
3083202ece2SStefano Zampini     if (!Bdense) {
3093202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
3103202ece2SStefano Zampini     } else {
3113202ece2SStefano Zampini       Bd = B;
3123202ece2SStefano Zampini     }
3133202ece2SStefano Zampini 
3143202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
3153202ece2SStefano Zampini       Mat fact;
3163202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
3173202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
3183202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
3193202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3203202ece2SStefano Zampini     } else {
32107b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
32207b1e237SStefano Zampini 
32307b1e237SStefano Zampini       if (ex) {
3243202ece2SStefano Zampini         Mat Ainvd;
3253202ece2SStefano Zampini 
3263202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3273202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3283202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
32907b1e237SStefano Zampini       } else {
33007b1e237SStefano Zampini         Vec         sol,rhs;
33107b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
33207b1e237SStefano Zampini         PetscInt    i,nrhs,n;
33307b1e237SStefano Zampini 
33407b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
33507b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
33607b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
33707b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
33807b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
33907b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
34007b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
34107b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
34207b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
34307b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
34407b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
34507b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
34607b1e237SStefano Zampini         }
34707b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
34807b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
34907b1e237SStefano Zampini       }
3503202ece2SStefano Zampini     }
3515ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3523202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3533202ece2SStefano Zampini     }
3545ec10c6aSStefano Zampini 
3555ec10c6aSStefano Zampini     if (!issym) {
3563202ece2SStefano Zampini       if (!Cdense) {
3573202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3583202ece2SStefano Zampini       } else {
3593202ece2SStefano Zampini         Cd = C;
3603202ece2SStefano Zampini       }
3615ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3623202ece2SStefano Zampini       if (!Cdense) {
3633202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3643202ece2SStefano Zampini       }
3655ec10c6aSStefano Zampini     } else {
3665ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3675ec10c6aSStefano Zampini       if (!Bdense) {
3685ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3695ec10c6aSStefano Zampini       }
3705ec10c6aSStefano Zampini     }
3715ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
372f11841e3SStefano Zampini   }
3733202ece2SStefano Zampini 
3743202ece2SStefano Zampini   if (D) {
3753202ece2SStefano Zampini     Mat       Dd;
3763202ece2SStefano Zampini     PetscBool Ddense;
3773202ece2SStefano Zampini 
3783202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3793202ece2SStefano Zampini     if (!Ddense) {
3803202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3813202ece2SStefano Zampini     } else {
3823202ece2SStefano Zampini       Dd = D;
3833202ece2SStefano Zampini     }
384f11841e3SStefano Zampini     if (n_I) {
3853202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
386f11841e3SStefano Zampini     } else {
387f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
388f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
389f11841e3SStefano Zampini       } else {
390f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
391f11841e3SStefano Zampini       }
392f11841e3SStefano Zampini     }
3933202ece2SStefano Zampini     if (!Ddense) {
3943202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3953202ece2SStefano Zampini     }
3963202ece2SStefano Zampini   } else {
3973202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3983202ece2SStefano Zampini   }
3993202ece2SStefano Zampini   PetscFunctionReturn(0);
4003202ece2SStefano Zampini }
40134a97f8cSStefano Zampini 
40291af6908SStefano 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)
403b1b3d7a2SStefano Zampini {
404be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
405be83ff47SStefano Zampini   Mat                    S_all;
4065a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
4075db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
408b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
409dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
410dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
411dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
4125a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
413683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
41408122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
41506a4b1faSStefano Zampini   PetscScalar            *Bwork;
4165a95e1ceSStefano Zampini   MPI_Comm               comm_n;
417f4f7d9d6SStefano Zampini   PetscBool              deluxe = PETSC_TRUE;
418f4f7d9d6SStefano Zampini   PetscBool              use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
4197f9db97bSStefano Zampini   PetscBool              print_schurs = PETSC_FALSE;
420b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
421b1b3d7a2SStefano Zampini 
422b1b3d7a2SStefano Zampini   PetscFunctionBegin;
423a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
424a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
425e62b6521Sstefano_zampini   /* convert matrix if needed */
426e62b6521Sstefano_zampini   if (Ain) {
427e62b6521Sstefano_zampini     PetscBool isseqaij;
428e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
429a64f4aa4SStefano Zampini     if (isseqaij) {
430a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
431a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4323301b35fSStefano Zampini     } else {
433a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
434a64f4aa4SStefano Zampini     }
435a64f4aa4SStefano Zampini   }
4363301b35fSStefano Zampini 
437a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
438a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
439df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
440df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
441a64f4aa4SStefano Zampini   }
442a64f4aa4SStefano Zampini 
4435a95e1ceSStefano Zampini   /* preliminary checks */
444af25d912SStefano 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");
4455a95e1ceSStefano Zampini 
44688113c35SStefano Zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
44788113c35SStefano Zampini 
4487f9db97bSStefano Zampini   /* debug prints */
4497f9db97bSStefano Zampini   if (sub_schurs->debug) {
4507f9db97bSStefano Zampini     PetscMPIInt size,rank;
4517f9db97bSStefano Zampini     PetscInt    nr,*print_schurs_ranks;
4527f9db97bSStefano Zampini     PetscBool   flg;
4537f9db97bSStefano Zampini 
4547f9db97bSStefano Zampini     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);CHKERRQ(ierr);
4557f9db97bSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4567f9db97bSStefano Zampini     nr   = size;
4577f9db97bSStefano Zampini     ierr = PetscMalloc1(nr,&print_schurs_ranks);CHKERRQ(ierr);
4587f9db97bSStefano Zampini     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
4597f9db97bSStefano Zampini     ierr = PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);CHKERRQ(ierr);
4607f9db97bSStefano Zampini     if (!flg) print_schurs = PETSC_TRUE;
4617f9db97bSStefano Zampini     else {
4627f9db97bSStefano Zampini       for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
4637f9db97bSStefano Zampini     }
4647f9db97bSStefano Zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
4657f9db97bSStefano Zampini     ierr = PetscFree(print_schurs_ranks);CHKERRQ(ierr);
4667f9db97bSStefano Zampini   }
4677f9db97bSStefano Zampini 
4687f9db97bSStefano Zampini 
4695a95e1ceSStefano Zampini   /* restrict work on active processes */
470991c41b4SStefano Zampini   if (sub_schurs->restrict_comm) {
471991c41b4SStefano Zampini     PetscSubcomm subcomm;
472991c41b4SStefano Zampini     PetscMPIInt  color,rank;
473991c41b4SStefano Zampini 
4745a95e1ceSStefano Zampini     color = 0;
4755a95e1ceSStefano Zampini     if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4765a95e1ceSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4775a95e1ceSStefano Zampini     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4785a95e1ceSStefano Zampini     ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4795a95e1ceSStefano Zampini     ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4805a95e1ceSStefano Zampini     ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4815a95e1ceSStefano Zampini     ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4825a95e1ceSStefano Zampini     if (!sub_schurs->n_subs) {
4835a95e1ceSStefano Zampini       ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4845a95e1ceSStefano Zampini       PetscFunctionReturn(0);
4855a95e1ceSStefano Zampini     }
486991c41b4SStefano Zampini   } else {
487991c41b4SStefano Zampini     ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr);
488991c41b4SStefano Zampini   }
4895a95e1ceSStefano Zampini 
490b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
491df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
492a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4933301b35fSStefano Zampini     PetscBool isseqsbaij;
494a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4953301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4963301b35fSStefano Zampini     if (isseqsbaij) {
497a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
498a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
499a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
500a64f4aa4SStefano Zampini     } else {
501a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
502a64f4aa4SStefano Zampini       A_BB = tA_BB;
503a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
504a64f4aa4SStefano Zampini       A_IB = tA_IB;
505a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
506a64f4aa4SStefano Zampini       A_BI = tA_BI;
507f11841e3SStefano Zampini     }
508a58a30b4SStefano Zampini   } else {
5095a95e1ceSStefano Zampini     A_II = NULL;
5105a95e1ceSStefano Zampini     A_IB = NULL;
5115a95e1ceSStefano Zampini     A_BI = NULL;
5125a95e1ceSStefano Zampini     A_BB = NULL;
513b1b3d7a2SStefano Zampini   }
5145a95e1ceSStefano Zampini   S_all = NULL;
515b1b3d7a2SStefano Zampini 
516b1b3d7a2SStefano Zampini   /* determine interior problems */
5173dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
5183dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
519b1b3d7a2SStefano Zampini     PetscBT                touched;
520b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
521b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
522b1b3d7a2SStefano Zampini 
5231633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
524b1b3d7a2SStefano Zampini     /* get sizes */
525b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
526b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
527b1b3d7a2SStefano Zampini 
528b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
529b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
530b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
531b1b3d7a2SStefano Zampini 
532b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
533b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
534b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
535b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
536b1b3d7a2SStefano Zampini     }
537b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
538b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
539b1b3d7a2SStefano Zampini 
540b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
541b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
542b1b3d7a2SStefano Zampini     n_prev_added = n_B;
543b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
544b1b3d7a2SStefano Zampini       PetscInt n_added;
545b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5466c4ed002SBarry 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);
547b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
548b1b3d7a2SStefano Zampini       n_prev_added = n_added;
549b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
550b1b3d7a2SStefano Zampini       if (!n_added) break;
551b1b3d7a2SStefano Zampini     }
552b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
553b1b3d7a2SStefano Zampini 
554883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
555a9b99552SStefano 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);
556b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
557a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
558883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
559df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
560b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
561b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
562a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
563b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
564b1b3d7a2SStefano Zampini 
565b1b3d7a2SStefano Zampini       /* II block */
5667dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
567b1b3d7a2SStefano Zampini     }
568b1b3d7a2SStefano Zampini   } else {
569b1b3d7a2SStefano Zampini     PetscInt n_I;
570b1b3d7a2SStefano Zampini 
571b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
572b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
573a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
574b1b3d7a2SStefano Zampini 
575b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
576df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
577b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
578b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
579b1b3d7a2SStefano Zampini 
580b1b3d7a2SStefano Zampini       /* II block is the same */
581b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
582b1b3d7a2SStefano Zampini       AE_II = A_II;
583b1b3d7a2SStefano Zampini     }
584b1b3d7a2SStefano Zampini   }
5855a95e1ceSStefano Zampini 
586883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5875a95e1ceSStefano Zampini   max_subset_size = 0;
588883469d8SStefano Zampini   local_size = 0;
5895a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5905a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5915a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
592883469d8SStefano Zampini     local_size += subset_size;
593883469d8SStefano Zampini   }
594883469d8SStefano Zampini 
595883469d8SStefano Zampini   /* Work arrays for local indices */
596883469d8SStefano Zampini   extra = 0;
597683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
598df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
599a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
600883469d8SStefano Zampini   }
601683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
602883469d8SStefano Zampini   if (extra) {
603883469d8SStefano Zampini     const PetscInt *idxs;
604a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
605883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
606a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
607883469d8SStefano Zampini   }
608883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
609dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
610dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
611883469d8SStefano Zampini 
612883469d8SStefano Zampini   /* Get local indices in local numbering */
613883469d8SStefano Zampini   local_size = 0;
6145a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
615883469d8SStefano Zampini     PetscInt j;
616883469d8SStefano Zampini     const    PetscInt *idxs;
617883469d8SStefano Zampini 
6185a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6195a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
620eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
621eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
622eb595f79SStefano Zampini     auxnum2[i] = subset_size;
623883469d8SStefano Zampini     /* subset indices in local numbering */
624883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
6255a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
626883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
627883469d8SStefano Zampini     local_size += subset_size;
628883469d8SStefano Zampini   }
629883469d8SStefano Zampini 
630f4f7d9d6SStefano Zampini   /* allocate extra workspace needed only for GETRI or SYTRF */
63111955456SStefano Zampini   use_potr = use_sytr = PETSC_FALSE;
63211955456SStefano Zampini   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
633f4f7d9d6SStefano Zampini     use_potr = PETSC_TRUE;
63411955456SStefano Zampini   } else if (sub_schurs->is_symmetric) {
63511955456SStefano Zampini     use_sytr = PETSC_TRUE;
63611955456SStefano Zampini   }
63711955456SStefano Zampini   if (local_size && !use_potr) {
63859ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
63959ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
640d2627357SStefano Zampini 
641d2627357SStefano Zampini     B_lwork = -1;
642d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
643d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
644f4f7d9d6SStefano Zampini     if (use_sytr) {
645f4f7d9d6SStefano Zampini       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
646f4f7d9d6SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
647f4f7d9d6SStefano Zampini     } else {
64859ac4de7SStefano Zampini       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
649d2627357SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
650f4f7d9d6SStefano Zampini     }
651f4f7d9d6SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
652d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
653d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
654056290a2SStefano Zampini   } else {
655056290a2SStefano Zampini     Bwork = NULL;
656056290a2SStefano Zampini     pivots = NULL;
657d2627357SStefano Zampini   }
658d2627357SStefano Zampini 
659d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
660dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
661dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
662dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
663dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6646583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
665dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
666dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
667dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6686c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
669dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
670e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
671d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
672d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
673d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
674d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6752972d61bSStefano Zampini 
6765a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6775a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6785a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6795a95e1ceSStefano Zampini 
6805a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6815a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
682*6080607fSStefano Zampini     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
6835a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
684b1b3d7a2SStefano Zampini   }
685b1b3d7a2SStefano Zampini 
68672b8c272SStefano Zampini   if (change) {
68772b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
68872b8c272SStefano Zampini     IS                     change_primal_B;
68972b8c272SStefano Zampini     IS                     change_primal_all;
69072b8c272SStefano Zampini 
691b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
692b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
693b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
69472b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
69572b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
69672b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
697b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
69872b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
69972b8c272SStefano Zampini     }
70072b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
70172b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
70272b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
70372b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
70472b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
705b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
70672b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
70772b8c272SStefano Zampini       Mat change_sub;
70872b8c272SStefano Zampini 
709b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
71072b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
71172b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
71272b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
7137dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
71472b8c272SStefano Zampini       } else {
71572b8c272SStefano Zampini         Mat change_subt;
7167dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
71772b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
71872b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
71972b8c272SStefano Zampini       }
72072b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
72172b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
722e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
723e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
72472b8c272SStefano Zampini     }
72572b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
72672b8c272SStefano Zampini   }
72772b8c272SStefano Zampini 
7285a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
7295a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
7305a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
7315a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
7325a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
7335a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
734aa83b6aeSStefano Zampini   }
735b1b3d7a2SStefano Zampini 
7365a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
737be83ff47SStefano Zampini   F = NULL;
738d943a642SStefano Zampini   if (!sub_schurs->schur_explicit) {
739d943a642SStefano Zampini     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
740d943a642SStefano Zampini        it is not efficient, unless the economic version of the scaling is used */
7415a95e1ceSStefano Zampini     Mat         S_Ej_expl;
7425a95e1ceSStefano Zampini     PetscScalar *work;
7435a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
7445a95e1ceSStefano Zampini     PetscBool   Sdense;
7455a95e1ceSStefano Zampini 
7465a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7475a95e1ceSStefano Zampini     local_size = 0;
748b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7495a95e1ceSStefano Zampini       IS  is_subset_B;
7505a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7515a95e1ceSStefano Zampini 
7525a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7535a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7545a95e1ceSStefano Zampini       /* EE block */
7557dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7565a95e1ceSStefano Zampini       /* IE block */
7577dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7585a95e1ceSStefano Zampini       /* EI block */
759d943a642SStefano Zampini       if (sub_schurs->is_symmetric) {
7605a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
761d943a642SStefano Zampini       } else if (sub_schurs->is_hermitian) {
762d943a642SStefano Zampini         ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7635a95e1ceSStefano Zampini       } else {
7647dae84e0SHong Zhang         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7655a95e1ceSStefano Zampini       }
766a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7675a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7685a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7695a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7705a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
771b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
772b1b3d7a2SStefano Zampini         KSP ksp;
773b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7745a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
775b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
776b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
777b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
778b1b3d7a2SStefano Zampini         KSPType   ksp_type;
779b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7805a95e1ceSStefano Zampini         PetscBool ispcnone;
781b1b3d7a2SStefano Zampini 
782b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7835a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
784b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
785b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
786b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
787b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7885a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7895a95e1ceSStefano Zampini         if (!ispcnone) {
7905a95e1ceSStefano Zampini           PCType pc_type;
791b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
792b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7935a95e1ceSStefano Zampini         } else {
7945a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7955a95e1ceSStefano Zampini         }
796b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
797b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
7983ca39a21SBarry Smith           MatSolverType solver = NULL;
799ea799195SBarry Smith           ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr);
800b1b3d7a2SStefano Zampini           if (solver) {
8013ca39a21SBarry Smith             ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr);
802b1b3d7a2SStefano Zampini           }
803b1b3d7a2SStefano Zampini         }
804b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
805b1b3d7a2SStefano Zampini       }
8065a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
8075a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
808d943a642SStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
8095a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
8105a95e1ceSStefano Zampini       if (Sdense) {
8115a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
8125a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
813b1b3d7a2SStefano Zampini         }
8145a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8156c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
8165a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
817a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
8185a95e1ceSStefano Zampini       local_size += subset_size;
8195a95e1ceSStefano Zampini     }
8205a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
821b1b3d7a2SStefano Zampini     /* free */
822b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
823b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
8245a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
825883469d8SStefano Zampini   } else {
8265cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
8279d54b7f4SStefano Zampini     Vec         Dall = NULL;
828ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
8295cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
8306dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
8313fc34f97SStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
8323fc34f97SStefano Zampini     PetscBool   schur_has_vertices,factor_workaround;
83311955456SStefano Zampini     PetscBool   use_cholesky;
834883469d8SStefano Zampini 
835683d3df6SStefano Zampini     /* get sizes */
83681ea8064SStefano Zampini     n_I = 0;
83781ea8064SStefano Zampini     if (is_I_layer) {
838a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
83981ea8064SStefano Zampini     }
840683d3df6SStefano Zampini     economic = PETSC_FALSE;
841683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
842683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
843683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
8449d54b7f4SStefano Zampini     size_active_schur = local_size;
8459d54b7f4SStefano Zampini 
846f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
8479d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
8489d54b7f4SStefano Zampini       const PetscScalar *array;
8499d54b7f4SStefano Zampini       PetscScalar       *array2;
8509d54b7f4SStefano Zampini       const PetscInt    *idxs;
8519d54b7f4SStefano Zampini       PetscInt          i;
8529d54b7f4SStefano Zampini 
8539d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8549d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8559d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8569d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8579d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8589d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8599d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8609d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8619d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8629d54b7f4SStefano Zampini     }
863d62866d3SStefano Zampini 
864683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
8653fc34f97SStefano Zampini     factor_workaround = PETSC_FALSE;
8663fc34f97SStefano Zampini     schur_has_vertices = PETSC_FALSE;
867683d3df6SStefano Zampini     cum = n_I+size_active_schur;
868683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
869683d3df6SStefano Zampini       const PetscInt* idxs;
870683d3df6SStefano Zampini       PetscInt        n_dir;
871683d3df6SStefano Zampini 
872683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
873683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
874683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
875683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
876683d3df6SStefano Zampini       cum += n_dir;
8773fc34f97SStefano Zampini       factor_workaround = PETSC_TRUE;
878d62866d3SStefano Zampini     }
879683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
880367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
881683d3df6SStefano Zampini       PetscInt n_v;
882683d3df6SStefano Zampini 
883683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
884683d3df6SStefano Zampini       if (n_v) {
885683d3df6SStefano Zampini         const PetscInt* idxs;
886683d3df6SStefano Zampini 
887683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
888683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
889683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
890683d3df6SStefano Zampini         cum += n_v;
891683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
8923fc34f97SStefano Zampini         schur_has_vertices = PETSC_TRUE;
893683d3df6SStefano Zampini       }
894683d3df6SStefano Zampini     }
895683d3df6SStefano Zampini     size_schur = cum - n_I;
896683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
897683d3df6SStefano Zampini     if (cum == n) {
898683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
899683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
900683d3df6SStefano Zampini     } else {
9017dae84e0SHong Zhang       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
902683d3df6SStefano Zampini     }
903e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
904e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
905ca92afb2SStefano Zampini 
906ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
907367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
908ca92afb2SStefano Zampini     if (benign_n) {
909ca92afb2SStefano Zampini       Vec                    v;
910ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
911ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
9125cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
9135cbda25cSStefano Zampini       PetscInt               sizeA;
914ca92afb2SStefano Zampini 
915ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
916ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
917ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
918ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
919ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
9205cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
9215cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
922282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
9235cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9245cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
925ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
926ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
927ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
9285cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
929ca92afb2SStefano Zampini         PetscScalar    *array;
930ca92afb2SStefano Zampini         const PetscInt *idxs;
931ca92afb2SStefano Zampini         PetscInt       j,nz;
932ca92afb2SStefano Zampini 
933ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
934ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
935ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9365cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
937ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9385cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9395cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
9405cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
941ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
94222db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
94322db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
94422db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
94522db5ddcSStefano Zampini #else
94622db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
94722db5ddcSStefano Zampini #endif
94822db5ddcSStefano Zampini         }
949ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
950ca92afb2SStefano Zampini       }
9515cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9525cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
953ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
954ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
955ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
956ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
957ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
958ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
959ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
960ca92afb2SStefano Zampini     }
96111955456SStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr);
9623301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9633301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
964883469d8SStefano Zampini 
96511955456SStefano Zampini     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
96611955456SStefano Zampini     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
96711955456SStefano Zampini 
968683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
969d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
970d47842beSStefano Zampini 
971f4f7d9d6SStefano Zampini     if (n_I) {
9720aa714b2SStefano Zampini       IS is_schur;
9735a05ddb0SStefano Zampini 
97411955456SStefano Zampini       if (use_cholesky) {
97588113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
976883469d8SStefano Zampini       } else {
97788113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
978883469d8SStefano Zampini       }
9798c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9808c09ecd8Sstefano_zampini 
981883469d8SStefano Zampini       /* subsets ordered last */
9825a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9835a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9845a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
985883469d8SStefano Zampini 
986883469d8SStefano Zampini       /* factorization step */
98711955456SStefano Zampini       if (use_cholesky) {
988883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
989be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
990be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
991be83ff47SStefano Zampini #endif
992883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
993a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
994883469d8SStefano Zampini       } else {
995883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
996be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
997be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
998be83ff47SStefano Zampini #endif
999883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1000883469d8SStefano Zampini       }
1001ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
1002883469d8SStefano Zampini 
10037f9db97bSStefano Zampini       if (print_schurs) {
100411955456SStefano Zampini         PetscViewer viewer;
100511955456SStefano Zampini         char        filename[256];
100611955456SStefano Zampini         Mat         S;
100711955456SStefano Zampini         IS          is;
100811955456SStefano Zampini 
10097f9db97bSStefano Zampini         ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
101011955456SStefano Zampini         ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
101111955456SStefano Zampini         ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
101211955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
101311955456SStefano Zampini         ierr = MatView(A,viewer);CHKERRQ(ierr);
101411955456SStefano Zampini         ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
101511955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr);
101611955456SStefano Zampini         ierr = MatView(S,viewer);CHKERRQ(ierr);
101711955456SStefano Zampini         ierr = MatDestroy(&S);CHKERRQ(ierr);
101811955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr);
101911955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr);
102011955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
102111955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
102211955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr);
102311955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr);
102411955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
102511955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
102611955456SStefano Zampini         ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
102711955456SStefano Zampini       }
102811955456SStefano Zampini 
1029883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
10307c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1031b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
1032b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
1033b3cb21ddSStefano Zampini 
1034d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
1035683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1036683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
103703dfb2d7SStefano Zampini       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
103803dfb2d7SStefano Zampini         reuse_solvers = factor_workaround = PETSC_FALSE;
103903dfb2d7SStefano Zampini 
1040df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
1041ca92afb2SStefano Zampini 
104272b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
1043ca92afb2SStefano Zampini       if (benign_n) {
10445cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1045ca92afb2SStefano Zampini         Vec         v;
10465cbda25cSStefano Zampini         PetscInt    sizeA;
1047ca92afb2SStefano Zampini 
1048ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
10495cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
10505cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1051ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1052ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1053ca92afb2SStefano Zampini #endif
1054ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1055ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1056ca92afb2SStefano Zampini #endif
10575cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10585cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1059ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
10605cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
106147484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
1062ca92afb2SStefano Zampini           const PetscInt *idxs;
1063ca92afb2SStefano Zampini           PetscInt       j,nz;
106447484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
1065ca92afb2SStefano Zampini 
10665cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
10675cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
10685cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1069ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1070ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1071282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
1072282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1073ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10745cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1075ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
107647484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
107747484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
107847484b83SStefano Zampini           B_k = 1;
107947484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
108047484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10815cbda25cSStefano Zampini           sum  = 1.;
108247484b83SStefano 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));
1083ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10845cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10855cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1086282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10875cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10885cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1089ca92afb2SStefano Zampini         }
1090a7414863SStefano Zampini         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1091a7414863SStefano Zampini           PetscInt k,j;
1092a7414863SStefano Zampini           for (k=0;k<size_schur;k++) {
1093a7414863SStefano Zampini             for (j=k;j<size_schur;j++) {
1094a7414863SStefano Zampini               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1095a7414863SStefano Zampini             }
1096a7414863SStefano Zampini           }
1097a7414863SStefano Zampini         }
1098a7414863SStefano Zampini 
10995cbda25cSStefano Zampini         /* restore defaults */
11005cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
11015cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
11025cbda25cSStefano Zampini #endif
11035cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
11045cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
11055cbda25cSStefano Zampini #endif
11065cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
11075cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1108ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1109ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1110ca92afb2SStefano Zampini       }
1111a3df083aSStefano Zampini       if (!reuse_solvers) {
1112a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1113a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1114a3df083aSStefano Zampini         }
1115a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
11165cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
11175cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1118a3df083aSStefano Zampini       }
1119df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1120be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1121683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1122166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1123df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1124be83ff47SStefano Zampini     }
1125be83ff47SStefano Zampini 
1126be83ff47SStefano Zampini     if (reuse_solvers) {
1127a00504b5SStefano Zampini       Mat                A_II,Afake;
112853892102SStefano Zampini       Vec                vec1_B;
1129df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
11303462e049SStefano Zampini       PetscInt           n_R;
1131d5574798SStefano Zampini 
1132df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1133df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1134e28d306cSStefano Zampini       } else {
1135df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1136d62866d3SStefano Zampini       }
1137df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1138be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1139d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1140d5574798SStefano Zampini       msolv_ctx->F = F;
1141683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1142683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1143683d3df6SStefano Zampini       {
1144683d3df6SStefano Zampini         PetscScalar *array;
1145683d3df6SStefano Zampini         PetscInt    n;
1146683d3df6SStefano Zampini 
1147683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1148683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1149683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1150683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1151683d3df6SStefano Zampini       }
11523fc34f97SStefano Zampini       msolv_ctx->has_vertices = schur_has_vertices;
1153d62866d3SStefano Zampini 
1154d62866d3SStefano Zampini       /* interior solver */
1155683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1156d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1157d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
115815579a77SStefano Zampini       ierr = PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");CHKERRQ(ierr);
1159d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
116015579a77SStefano Zampini       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1161df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1162df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1163d62866d3SStefano Zampini 
1164d62866d3SStefano Zampini       /* correction solver */
11653462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1166d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
116715579a77SStefano Zampini       ierr = PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");CHKERRQ(ierr);
1168d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
116915579a77SStefano Zampini       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1170df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1171df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
117253892102SStefano Zampini 
117353892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
117453892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
117553892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
11763fc34f97SStefano Zampini       if (!schur_has_vertices) {
117753892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
117853892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
117953892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
118053892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1181683d3df6SStefano Zampini       } else {
1182683d3df6SStefano Zampini         IS              is_B_all;
1183683d3df6SStefano Zampini         const PetscInt* idxs;
1184683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1185683d3df6SStefano Zampini 
1186683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1187683d3df6SStefano Zampini         dual = size_schur - n_v;
1188683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1189683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1190683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1191683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1192683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1193683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1194683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1195683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1196683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1197683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1198683d3df6SStefano Zampini       }
11993462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
12003462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
12018d692618SStefano Zampini       ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12028d692618SStefano Zampini       ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12033462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
12043462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1205683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1206ca92afb2SStefano Zampini 
1207ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1208ca92afb2SStefano Zampini       if (benign_n) {
12095cbda25cSStefano Zampini         PetscScalar *array;
12105cbda25cSStefano Zampini 
1211ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1212ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1213ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
12145cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1215282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
12165cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
12175cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
12185cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
12195cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1220ca92afb2SStefano Zampini       }
1221ada6e2d7SStefano Zampini     } else {
1222ada6e2d7SStefano Zampini       if (sub_schurs->reuse_solver) {
1223ada6e2d7SStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1224ada6e2d7SStefano Zampini       }
1225ada6e2d7SStefano Zampini       ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1226d5574798SStefano Zampini     }
122708122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
122853892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
12295db18549SStefano Zampini 
1230be83ff47SStefano Zampini     /* Work arrays */
1231be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
12325a94c880SStefano Zampini 
123304270d10SStefano Zampini     /* matrices for deluxe scaling and adaptive selection */
123412d906b1SStefano Zampini     if (compute_Stilda) {
12355a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
1236b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
12375a94c880SStefano Zampini       }
12389d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1239b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
12405a94c880SStefano Zampini       }
124112d906b1SStefano Zampini     }
1242d2627357SStefano Zampini 
1243be83ff47SStefano Zampini     /* S_Ej_all */
1244be83ff47SStefano Zampini     cum = cum2 = 0;
1245be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
124665d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
124765d8bf0aSStefano Zampini       PetscInt j;
124865d8bf0aSStefano Zampini 
12495a95e1ceSStefano Zampini       /* get S_E */
1250b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1251683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1252be83ff47SStefano Zampini         PetscInt k;
1253be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1254be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1255be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1256683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1257be83ff47SStefano Zampini           }
1258be83ff47SStefano Zampini         }
125906a4e24aSStefano Zampini       } else { /* just copy to workspace */
1260be83ff47SStefano Zampini         PetscInt k;
1261be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1262be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1263be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1264be83ff47SStefano Zampini           }
1265be83ff47SStefano Zampini         }
12669087bf02SStefano Zampini       }
12675a95e1ceSStefano Zampini       /* insert S_E values */
1268be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1269b7ab4a40SStefano Zampini       if (sub_schurs->change) {
12708760537fSStefano Zampini         Mat change_sub,SEj,T;
127172b8c272SStefano Zampini 
127272b8c272SStefano Zampini         /* change basis */
127372b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
127472b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12758760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
12768760537fSStefano Zampini           Mat T2;
12778760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
12788760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
12798760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
128072b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
12818760537fSStefano Zampini         } else {
12828760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
128372b8c272SStefano Zampini         }
128472b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
128572b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1286afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
128772b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
128872b8c272SStefano Zampini       }
12899d54b7f4SStefano Zampini       if (deluxe) {
12905a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1291683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1292862806e4SStefano Zampini         if (compute_Stilda) {
1293f4f7d9d6SStefano Zampini           PetscInt k;
1294f4f7d9d6SStefano Zampini 
129508122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
129608122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1297f4f7d9d6SStefano Zampini           if (use_potr) {
12985a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12992972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
13005a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
13012972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
13026c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
13036c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
13046c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
13056c3e6151SStefano Zampini               }
13066c3e6151SStefano Zampini             }
1307f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1308f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1309f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1310f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1311f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1312f4f7d9d6SStefano Zampini             for (k=0;k<subset_size;k++) {
1313f4f7d9d6SStefano Zampini               for (j=k;j<subset_size;j++) {
1314f4f7d9d6SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
1315f4f7d9d6SStefano Zampini               }
1316f4f7d9d6SStefano Zampini             }
1317d6462365SStefano Zampini           } else {
1318d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&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,work,&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);
13222972d61bSStefano Zampini           }
132343371fb9SStefano Zampini           ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr);
132408122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
13255a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13269087bf02SStefano Zampini         }
13279d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
13289d54b7f4SStefano Zampini         Mat         SEj;
13299d54b7f4SStefano Zampini         Vec         D;
13309d54b7f4SStefano Zampini         PetscScalar *array;
13319d54b7f4SStefano Zampini 
13329d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
13339d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
13349d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
13359d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1336f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
13379d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
13389d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
13399d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
13409d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13419d54b7f4SStefano Zampini       }
1342be83ff47SStefano Zampini       cum += subset_size;
1343be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1344be83ff47SStefano Zampini     }
1345be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
13464a6c6b0dSStefano Zampini 
1347df4d28bfSStefano Zampini     if (solver_S) {
13487c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
13494a6c6b0dSStefano Zampini     }
1350683d3df6SStefano Zampini 
1351683d3df6SStefano Zampini     schur_factor = NULL;
135245951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1353683d3df6SStefano Zampini 
13549d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
13554a6c6b0dSStefano Zampini         PetscInt j;
13564a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
13575a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13584a6c6b0dSStefano Zampini       } else {
1359683d3df6SStefano Zampini         Mat S_all_inv=NULL;
13603fc34f97SStefano Zampini         if (solver_S) {
1361683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1362683d3df6SStefano 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 */
13633fc34f97SStefano Zampini           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1364683d3df6SStefano Zampini             PetscScalar *data;
1365683d3df6SStefano Zampini             PetscInt     nd = 0;
13666dba178dSStefano Zampini 
1367f4f7d9d6SStefano Zampini             if (!use_potr) {
1368683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1369683d3df6SStefano Zampini             }
13707c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1371683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1372683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1373683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1374683d3df6SStefano Zampini             }
13753fc34f97SStefano Zampini 
13763fc34f97SStefano Zampini             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
13773fc34f97SStefano Zampini             if (schur_has_vertices) {
13783fc34f97SStefano Zampini               Mat          M;
13793fc34f97SStefano Zampini               PetscScalar *tdata;
13803fc34f97SStefano Zampini               PetscInt     nv = 0, news;
13813fc34f97SStefano Zampini 
13823fc34f97SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
13833fc34f97SStefano Zampini               news = size_active_schur + nv;
13843fc34f97SStefano Zampini               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1385683d3df6SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13863fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13873fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr);
13883fc34f97SStefano Zampini               }
13893fc34f97SStefano Zampini               for (i=0;i<nv;i++) {
13903fc34f97SStefano Zampini                 PetscInt k = i+size_active_schur;
13913fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13923fc34f97SStefano Zampini               }
13933fc34f97SStefano Zampini 
13943fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
13953fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13963fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13973fc34f97SStefano Zampini               /* save the factors */
13983fc34f97SStefano Zampini               cum = 0;
13993fc34f97SStefano Zampini               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
14003fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
14013fc34f97SStefano Zampini                 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1402683d3df6SStefano Zampini                 cum += size_active_schur - i;
1403683d3df6SStefano Zampini               }
14043fc34f97SStefano Zampini               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
14053fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
14063fc34f97SStefano Zampini               /* move back just the active dofs to the Schur complement */
14073fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
14083fc34f97SStefano Zampini                 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr);
14093fc34f97SStefano Zampini               }
14103fc34f97SStefano Zampini               ierr = PetscFree(tdata);CHKERRQ(ierr);
14113fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
14123fc34f97SStefano Zampini             } else { /* we can factorize and invert just the activedofs */
14133fc34f97SStefano Zampini               Mat         M;
14145002105bSStefano Zampini               PetscScalar *aux;
14153fc34f97SStefano Zampini 
14165002105bSStefano Zampini               ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr);
14175002105bSStefano Zampini               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
14183fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
14193fc34f97SStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
14203fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
14213fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
14223fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
14233fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
14245002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr);
14255002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
14265002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
14275002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr);
14285002105bSStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
14295002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
14305002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
14315002105bSStefano Zampini               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
14325002105bSStefano Zampini               ierr = PetscFree(aux);CHKERRQ(ierr);
14333fc34f97SStefano Zampini             }
1434683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
14353fc34f97SStefano Zampini           } else { /* use MatFactor calls to invert S */
14366dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
14377c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1438683d3df6SStefano Zampini           }
1439df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1440683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1441683d3df6SStefano Zampini           S_all_inv = S_all;
1442683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1443be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1444be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1445f4f7d9d6SStefano Zampini           if (use_potr) {
1446be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1447be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1448be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1449be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1450f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1451f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1452f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1453f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1454f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1455d6462365SStefano Zampini           } else {
1456d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1457d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1458d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1459d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1460be83ff47SStefano Zampini           }
146143371fb9SStefano Zampini           ierr = PetscLogFlops(1.0*size_schur*size_schur*size_schur);CHKERRQ(ierr);
1462be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1463683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1464be83ff47SStefano Zampini         }
1465be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1466be83ff47SStefano Zampini         cum = cum2 = 0;
1467683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1468be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1469be83ff47SStefano Zampini           PetscInt j;
1470862806e4SStefano Zampini 
1471862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1472be83ff47SStefano Zampini           /* get (St^-1)_E */
147372b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
147406a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1475a0b0af32SStefano Zampini           if (S_lower_triangular) {
1476be83ff47SStefano Zampini             PetscInt k;
1477b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1478be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1479be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1480be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
14816c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1482be83ff47SStefano Zampini                 }
1483be83ff47SStefano Zampini               }
148472b8c272SStefano Zampini             } else {
148572b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
148672b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
148772b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
148872b8c272SStefano Zampini                 }
148972b8c272SStefano Zampini               }
149072b8c272SStefano Zampini             }
149172b8c272SStefano Zampini           } else {
1492be83ff47SStefano Zampini             PetscInt k;
1493be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1494be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1495be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1496be83ff47SStefano Zampini               }
1497be83ff47SStefano Zampini             }
1498be83ff47SStefano Zampini           }
1499b7ab4a40SStefano Zampini           if (sub_schurs->change) {
15008760537fSStefano Zampini             Mat change_sub,SEj,T;
150172b8c272SStefano Zampini 
150272b8c272SStefano Zampini             /* change basis */
150372b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
150472b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
15058760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
15068760537fSStefano Zampini               Mat T2;
15078760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
15088760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
150972b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
15108760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
15118760537fSStefano Zampini             } else {
15128760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
151372b8c272SStefano Zampini             }
15148760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
15158760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
151672b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1517afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
151872b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
151972b8c272SStefano Zampini           }
1520be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
15215a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1522be83ff47SStefano Zampini           cum += subset_size;
1523be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1524883469d8SStefano Zampini         }
1525683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1526df4d28bfSStefano Zampini         if (solver_S) {
15273fc34f97SStefano Zampini           if (schur_has_vertices) {
15283fc34f97SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
15293fc34f97SStefano Zampini           } else {
15307c2f51b8SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
15315db18549SStefano Zampini           }
15323fc34f97SStefano Zampini         }
1533683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1534683d3df6SStefano Zampini       }
1535683d3df6SStefano Zampini 
15363fc34f97SStefano Zampini       /* move back factors if needed */
15373fc34f97SStefano Zampini       if (schur_has_vertices) {
1538683d3df6SStefano Zampini         Mat      S_tmp;
15393fc34f97SStefano Zampini         PetscInt nd = 0;
1540683d3df6SStefano Zampini 
15416542c05fSStefano Zampini         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
15427c2f51b8SStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1543f4f7d9d6SStefano Zampini         if (use_potr) {
1544683d3df6SStefano Zampini           PetscScalar *data;
1545683d3df6SStefano Zampini 
1546683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
15473fc34f97SStefano Zampini           ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1548683d3df6SStefano Zampini 
1549683d3df6SStefano Zampini           if (S_lower_triangular) {
1550683d3df6SStefano Zampini             cum = 0;
1551683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1552683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1553683d3df6SStefano Zampini               cum += size_active_schur-i;
1554683d3df6SStefano Zampini 	    }
1555683d3df6SStefano Zampini           } else {
1556683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1557683d3df6SStefano Zampini           }
1558683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1559683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1560683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1561683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1562683d3df6SStefano Zampini 	    }
1563683d3df6SStefano Zampini           }
15646dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1565683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1566683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
15676c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1568683d3df6SStefano Zampini           }
1569683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
15706542c05fSStefano Zampini         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
15713fc34f97SStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
15729087bf02SStefano Zampini       }
1573367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1574367aa537SStefano Zampini       PetscScalar *data;
1575367aa537SStefano Zampini       PetscInt    nd = 0;
1576367aa537SStefano Zampini 
1577367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1578367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1579367aa537SStefano Zampini       }
15807c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1581367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1582367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1583367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1584367aa537SStefano Zampini       }
1585367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1586367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
15876c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1588367aa537SStefano Zampini       }
1589367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
15903fc34f97SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
15914a6c6b0dSStefano Zampini     }
15924a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1593683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
15949d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
15954a6c6b0dSStefano Zampini   }
15965a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1597be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1598a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1599a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1600a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1601a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1602a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
16036afe12f5SStefano Zampini 
16046afe12f5SStefano Zampini   /* speed up matrix assembly */
16056afe12f5SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&nnz);CHKERRQ(ierr);
16066afe12f5SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
16076afe12f5SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);CHKERRQ(ierr);
16086afe12f5SStefano Zampini   }
16096afe12f5SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);CHKERRQ(ierr);
16106afe12f5SStefano Zampini   ierr = MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
16115db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16125db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16135a95e1ceSStefano Zampini   if (compute_Stilda) {
16146afe12f5SStefano Zampini     ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
16155a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16165a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16179d54b7f4SStefano Zampini     if (deluxe) {
16186afe12f5SStefano Zampini       ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
16195a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16205a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162108122e43SStefano Zampini     }
16229d54b7f4SStefano Zampini   }
16236afe12f5SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
16246afe12f5SStefano Zampini 
16255db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
16265db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1627487b449aSStefano Zampini   ierr = MatAssemblyBegin(work_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1628487b449aSStefano Zampini   ierr = MatAssemblyEnd(work_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16293927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1630487b449aSStefano Zampini   ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
16315a95e1ceSStefano Zampini 
16325db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
16337dae84e0SHong Zhang   ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
163441fd5443SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
163541fd5443SStefano Zampini     ierr = MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
163641fd5443SStefano Zampini   } else {
163704270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
163841fd5443SStefano Zampini   }
163908122e43SStefano Zampini 
1640f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
16415a95e1ceSStefano Zampini   if (compute_Stilda) {
16425a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1643487b449aSStefano Zampini     ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
164404270d10SStefano Zampini     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
164504270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
16469d54b7f4SStefano Zampini     if (deluxe) {
16475a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1648487b449aSStefano Zampini       ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
164904270d10SStefano Zampini       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
165004270d10SStefano Zampini       ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
16519d54b7f4SStefano Zampini     } else {
16529d54b7f4SStefano Zampini       PetscScalar *array;
16539d54b7f4SStefano Zampini       PetscInt    cum;
16549d54b7f4SStefano Zampini 
16559d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
16569d54b7f4SStefano Zampini       cum = 0;
16579d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
16589d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
16599d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
16609d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1661f4f7d9d6SStefano Zampini         if (use_potr) {
16629d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
16639d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
16649d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
16659d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1666f4f7d9d6SStefano Zampini         } else if (use_sytr) {
1667f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1668f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1669f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1670f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1671f4f7d9d6SStefano Zampini         } else {
1672f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1673f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1674f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1675f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1676f4f7d9d6SStefano Zampini         }
167743371fb9SStefano Zampini         ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr);
16789d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
16799d54b7f4SStefano Zampini         cum += subset_size*subset_size;
16809d54b7f4SStefano Zampini       }
16819d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
16829d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
16839d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
16849d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
16859d54b7f4SStefano Zampini     }
168608122e43SStefano Zampini   }
168704270d10SStefano Zampini   ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr);
16887f9db97bSStefano Zampini   if (print_schurs) {
168911955456SStefano Zampini     PetscViewer viewer;
169011955456SStefano Zampini     char        filename[256];
169111955456SStefano Zampini     PetscInt    cum;
169211955456SStefano Zampini 
16937f9db97bSStefano Zampini     ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_mats_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
169411955456SStefano Zampini     ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
169511955456SStefano Zampini     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
169611955456SStefano Zampini     if (sub_schurs->S_Ej_all) {
169711955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr);
169811955456SStefano Zampini       ierr = MatView(sub_schurs->S_Ej_all,viewer);CHKERRQ(ierr);
169911955456SStefano Zampini     }
170011955456SStefano Zampini     if (sub_schurs->sum_S_Ej_all) {
170111955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr);
170211955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_all,viewer);CHKERRQ(ierr);
170311955456SStefano Zampini     }
170411955456SStefano Zampini     if (sub_schurs->sum_S_Ej_inv_all) {
170511955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr);
170611955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_inv_all,viewer);CHKERRQ(ierr);
170711955456SStefano Zampini     }
170811955456SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
170911955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr);
171011955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);CHKERRQ(ierr);
171111955456SStefano Zampini     }
171211955456SStefano Zampini     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
171311955456SStefano Zampini       IS   is;
171411955456SStefano Zampini       char name[16];
171511955456SStefano Zampini 
17167f9db97bSStefano Zampini       ierr = PetscSNPrintf(name,sizeof(name),"IE%D",i);CHKERRQ(ierr);
171711955456SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
171811955456SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr);
171911955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr);
172011955456SStefano Zampini       ierr = ISView(is,viewer);CHKERRQ(ierr);
172111955456SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
172211955456SStefano Zampini       cum += subset_size;
172311955456SStefano Zampini     }
172411955456SStefano Zampini     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
172511955456SStefano Zampini   }
17263202ece2SStefano Zampini 
17275a95e1ceSStefano Zampini   /* free workspace */
17285a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
172906a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1730a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
17313202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1732dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
17335a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1734b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1735b1b3d7a2SStefano Zampini }
1736b1b3d7a2SStefano Zampini 
173788113c35SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1738b1b3d7a2SStefano Zampini {
17399bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
17405a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
174188113c35SStefano Zampini   PetscBool       is_sorted,ispetsc;
1742b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1743b1b3d7a2SStefano Zampini 
1744b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1745b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
17466c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1747b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
17486c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1749b1b3d7a2SStefano Zampini 
1750b1b3d7a2SStefano Zampini   /* reset any previous data */
1751b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1752b1b3d7a2SStefano Zampini 
17535a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
17549bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1755b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
175608122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1757b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1758b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
17598b6046baSStefano Zampini     if (copycc) {
17608b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
17618b6046baSStefano Zampini     } else {
1762c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1763b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1764b1b3d7a2SStefano Zampini     }
17658b6046baSStefano Zampini   }
1766b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
17678b6046baSStefano Zampini     if (copycc) {
17688b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
17698b6046baSStefano Zampini     } else {
1770c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1771b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
17728b6046baSStefano Zampini     }
177308122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1774b1b3d7a2SStefano Zampini   }
1775c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1776c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1777c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1778d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1779d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1780b1b3d7a2SStefano Zampini 
1781df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
178288113c35SStefano Zampini   ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr);
1783883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
178488113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);CHKERRQ(ierr);
178588113c35SStefano Zampini #elif defined(PETSC_HAVE_MKL_PARDISO)
178688113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);CHKERRQ(ierr);
178788113c35SStefano Zampini #else
178888113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);CHKERRQ(ierr);
1789df4d28bfSStefano Zampini #endif
179088113c35SStefano Zampini #if defined(PETSC_USE_COMPLEX)
179188113c35SStefano Zampini   sub_schurs->is_hermitian  = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
179288113c35SStefano Zampini #else
179388113c35SStefano Zampini   sub_schurs->is_hermitian  = PETSC_TRUE;
1794883469d8SStefano Zampini #endif
179588113c35SStefano Zampini   sub_schurs->is_posdef     = PETSC_TRUE;
179611955456SStefano Zampini   sub_schurs->is_symmetric  = PETSC_TRUE;
17977f9db97bSStefano Zampini   sub_schurs->debug         = PETSC_FALSE;
1798991c41b4SStefano Zampini   sub_schurs->restrict_comm = PETSC_FALSE;
179988113c35SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
180088113c35SStefano Zampini   ierr = PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);CHKERRQ(ierr);
180111955456SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr);
180288113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
180388113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
1804991c41b4SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);CHKERRQ(ierr);
18057f9db97bSStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);CHKERRQ(ierr);
180688113c35SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
180788113c35SStefano Zampini   ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);CHKERRQ(ierr);
1808168c3995SStefano Zampini   sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1809b1b3d7a2SStefano Zampini 
181011955456SStefano Zampini   /* for reals, symmetric and hermitian are synonims */
181111955456SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
181211955456SStefano Zampini   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
181311955456SStefano Zampini   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
181411955456SStefano Zampini #endif
181511955456SStefano Zampini 
1816b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1817b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1818b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1819b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
18205db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
18215db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
18225db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
18235db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
18245a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1825b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1826b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1827b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
182808122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1829b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1830b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1831b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1832b1b3d7a2SStefano Zampini }
1833b1b3d7a2SStefano Zampini 
183434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
183534a97f8cSStefano Zampini {
183634a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
183734a97f8cSStefano Zampini   PetscErrorCode  ierr;
183834a97f8cSStefano Zampini 
183934a97f8cSStefano Zampini   PetscFunctionBegin;
184034a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
18415ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
184234a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
184334a97f8cSStefano Zampini   PetscFunctionReturn(0);
184434a97f8cSStefano Zampini }
184534a97f8cSStefano Zampini 
184634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
184734a97f8cSStefano Zampini {
184834a97f8cSStefano Zampini   PetscInt       i;
184934a97f8cSStefano Zampini   PetscErrorCode ierr;
185034a97f8cSStefano Zampini 
185134a97f8cSStefano Zampini   PetscFunctionBegin;
1852aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
185388113c35SStefano Zampini   ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr);
18541e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1855b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1856b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1857b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
18585db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
18595db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
186041c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
186141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
186208122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1863a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
18645db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1865d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1866d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
186708122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
186834a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1869b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
187034a97f8cSStefano Zampini   }
18715ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1872b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
18733dc780c3SStefano Zampini   }
1874df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1875df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1876d62866d3SStefano Zampini   }
1877df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
187872b8c272SStefano Zampini   if (sub_schurs->change) {
187972b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
188072b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1881b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
188272b8c272SStefano Zampini     }
188372b8c272SStefano Zampini   }
188472b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1885b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
188634a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
188734a97f8cSStefano Zampini   PetscFunctionReturn(0);
188834a97f8cSStefano Zampini }
188934a97f8cSStefano Zampini 
1890aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1891aea80f77Sstefano_zampini {
1892aea80f77Sstefano_zampini   PetscErrorCode ierr;
1893aea80f77Sstefano_zampini 
1894aea80f77Sstefano_zampini   PetscFunctionBegin;
1895aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1896aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1897aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1898aea80f77Sstefano_zampini }
1899aea80f77Sstefano_zampini 
19002a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
190134a97f8cSStefano Zampini {
190234a97f8cSStefano Zampini   PetscInt       i,j,n;
190334a97f8cSStefano Zampini   PetscErrorCode ierr;
190434a97f8cSStefano Zampini 
190534a97f8cSStefano Zampini   PetscFunctionBegin;
190634a97f8cSStefano Zampini   n = 0;
190734a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
190834a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
190934a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
191034a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
191134a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
191234a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
191334a97f8cSStefano Zampini         queue_tip[n] = dof;
191434a97f8cSStefano Zampini         n++;
191534a97f8cSStefano Zampini       }
191634a97f8cSStefano Zampini     }
191734a97f8cSStefano Zampini   }
191834a97f8cSStefano Zampini   *n_added = n;
191934a97f8cSStefano Zampini   PetscFunctionReturn(0);
192034a97f8cSStefano Zampini }
1921