xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision cc17cd1eb2c95d3bd5cf11ca769fdbee5b9ef269)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 #include <../src/mat/impls/dense/seq/dense.h>
4 #include <petscblaslapack.h>
5 
6 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
7 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
8 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
9 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
10 
11 /* if v2 is not present, correction is done in-place */
12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13 {
14   PetscScalar    *array;
15   PetscScalar    *array2;
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   if (!ctx->benign_n) PetscFunctionReturn(0);
20   if (sol && full) {
21     PetscInt n_I,size_schur;
22 
23     /* get sizes */
24     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
25     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
26     n_I = n_I - size_schur;
27     /* get schur sol from array */
28     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
29     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
30     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
31     /* apply interior sol correction */
32     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
33     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
34     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
35   }
36   if (v2) {
37     PetscInt nl;
38 
39     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
40     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
41     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
42     ierr = PetscArraycpy(array2,array,nl);CHKERRQ(ierr);
43   } else {
44     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
45     array2 = array;
46   }
47   if (!sol) { /* change rhs */
48     PetscInt n;
49     for (n=0;n<ctx->benign_n;n++) {
50       PetscScalar    sum = 0.;
51       const PetscInt *cols;
52       PetscInt       nz,i;
53 
54       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
55       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
56       for (i=0;i<nz-1;i++) sum += array[cols[i]];
57 #if defined(PETSC_USE_COMPLEX)
58       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
59 #else
60       sum = -sum/nz;
61 #endif
62       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
63       ctx->benign_save_vals[n] = array2[cols[nz-1]];
64       array2[cols[nz-1]] = sum;
65       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
66     }
67   } else {
68     PetscInt n;
69     for (n=0;n<ctx->benign_n;n++) {
70       PetscScalar    sum = 0.;
71       const PetscInt *cols;
72       PetscInt       nz,i;
73       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
74       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
75       for (i=0;i<nz-1;i++) sum += array[cols[i]];
76 #if defined(PETSC_USE_COMPLEX)
77       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
78 #else
79       sum = -sum/nz;
80 #endif
81       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
82       array2[cols[nz-1]] = ctx->benign_save_vals[n];
83       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
84     }
85   }
86   if (v2) {
87     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
88     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
89   } else {
90     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
91   }
92   if (!sol && full) {
93     Vec      usedv;
94     PetscInt n_I,size_schur;
95 
96     /* get sizes */
97     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
98     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
99     n_I = n_I - size_schur;
100     /* compute schur rhs correction */
101     if (v2) {
102       usedv = v2;
103     } else {
104       usedv = v;
105     }
106     /* apply schur rhs correction */
107     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
108     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
109     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
110     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
111     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
112     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
113   }
114   PetscFunctionReturn(0);
115 }
116 
117 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118 {
119   PCBDDCReuseSolvers ctx;
120   PetscBool          copy = PETSC_FALSE;
121   PetscErrorCode     ierr;
122 
123   PetscFunctionBegin;
124   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
125   if (full) {
126 #if defined(PETSC_HAVE_MUMPS)
127     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
128 #endif
129 #if defined(PETSC_HAVE_MKL_PARDISO)
130     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
131 #endif
132     copy = ctx->has_vertices;
133   } else { /* interior solver */
134 #if defined(PETSC_HAVE_MUMPS)
135     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
136 #endif
137 #if defined(PETSC_HAVE_MKL_PARDISO)
138     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
139 #endif
140     copy = PETSC_TRUE;
141   }
142   /* copy rhs into factored matrix workspace */
143   if (copy) {
144     PetscInt    n;
145     PetscScalar *array,*array_solver;
146 
147     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
148     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
149     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
150     ierr = PetscArraycpy(array_solver,array,n);CHKERRQ(ierr);
151     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
152     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
153 
154     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
155     if (transpose) {
156       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
157     } else {
158       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
159     }
160     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
161 
162     /* get back data to caller worskpace */
163     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
164     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
165     ierr = PetscArraycpy(array,array_solver,n);CHKERRQ(ierr);
166     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
167     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
168   } else {
169     if (ctx->benign_n) {
170       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
171       if (transpose) {
172         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
173       } else {
174         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
175       }
176       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
177     } else {
178       if (transpose) {
179         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
180       } else {
181         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
182       }
183     }
184   }
185   /* restore defaults */
186 #if defined(PETSC_HAVE_MUMPS)
187   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
188 #endif
189 #if defined(PETSC_HAVE_MKL_PARDISO)
190   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
191 #endif
192   PetscFunctionReturn(0);
193 }
194 
195 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196 {
197   PetscErrorCode   ierr;
198 
199   PetscFunctionBegin;
200   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
201   PetscFunctionReturn(0);
202 }
203 
204 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205 {
206   PetscErrorCode   ierr;
207 
208   PetscFunctionBegin;
209   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214 {
215   PetscErrorCode   ierr;
216 
217   PetscFunctionBegin;
218   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223 {
224   PetscErrorCode   ierr;
225 
226   PetscFunctionBegin;
227   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
228   PetscFunctionReturn(0);
229 }
230 
231 static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer)
232 {
233   PCBDDCReuseSolvers ctx;
234   PetscBool          iascii;
235   PetscErrorCode     ierr;
236 
237   PetscFunctionBegin;
238   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
239   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
240   if (iascii) {
241     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
242   }
243   ierr = MatView(ctx->F,viewer);CHKERRQ(ierr);
244   if (iascii) {
245     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
246   }
247   PetscFunctionReturn(0);
248 }
249 
250 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
251 {
252   PetscInt       i;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
257   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
258   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
259   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
260   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
261   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
262   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
263   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
264   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
265   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
266   for (i=0;i<reuse->benign_n;i++) {
267     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
268   }
269   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
270   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
271   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
272   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
273   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
274   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
279 {
280   Mat            B, C, D, Bd, Cd, AinvBd;
281   KSP            ksp;
282   PC             pc;
283   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
284   PetscReal      fill = 2.0;
285   PetscInt       n_I;
286   PetscMPIInt    size;
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
291   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
292   if (reuse == MAT_REUSE_MATRIX) {
293     PetscBool Sdense;
294 
295     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
296     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
297   }
298   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
299   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
300   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
301   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
302   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
303   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
304   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
305   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
306   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
307   if (n_I) {
308     if (!Bdense) {
309       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
310     } else {
311       Bd = B;
312     }
313 
314     if (isLU || isILU || isCHOL) {
315       Mat fact;
316       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
317       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
318       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
319       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
320     } else {
321       PetscBool ex = PETSC_TRUE;
322 
323       if (ex) {
324         Mat Ainvd;
325 
326         ierr = PCComputeOperator(pc, MATDENSE, &Ainvd);CHKERRQ(ierr);
327         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
328         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
329       } else {
330         Vec         sol,rhs;
331         PetscScalar *arrayrhs,*arraysol;
332         PetscInt    i,nrhs,n;
333 
334         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
335         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
336         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
337         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
338         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
339         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
340         for (i=0;i<nrhs;i++) {
341           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
342           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
343           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
344           ierr = VecResetArray(rhs);CHKERRQ(ierr);
345           ierr = VecResetArray(sol);CHKERRQ(ierr);
346         }
347         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
348         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
349       }
350     }
351     if (!Bdense & !issym) {
352       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
353     }
354 
355     if (!issym) {
356       if (!Cdense) {
357         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
358       } else {
359         Cd = C;
360       }
361       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
362       if (!Cdense) {
363         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
364       }
365     } else {
366       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
367       if (!Bdense) {
368         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
369       }
370     }
371     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
372   }
373 
374   if (D) {
375     Mat       Dd;
376     PetscBool Ddense;
377 
378     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
379     if (!Ddense) {
380       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
381     } else {
382       Dd = D;
383     }
384     if (n_I) {
385       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
386     } else {
387       if (reuse == MAT_INITIAL_MATRIX) {
388         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
389       } else {
390         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
391       }
392     }
393     if (!Ddense) {
394       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
395     }
396   } else {
397     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 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)
403 {
404   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
405   Mat                    S_all;
406   Vec                    gstash,lstash;
407   VecScatter             sstash;
408   IS                     is_I,is_I_layer;
409   IS                     all_subsets,all_subsets_mult,all_subsets_n;
410   PetscScalar            *stasharray,*Bwork;
411   PetscInt               *nnz,*all_local_idx_N;
412   PetscInt               *auxnum1,*auxnum2;
413   PetscInt               i,subset_size,max_subset_size;
414   PetscInt               n_B,extra,local_size,global_size;
415   PetscInt               local_stash_size;
416   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
417   MPI_Comm               comm_n;
418   PetscBool              deluxe = PETSC_TRUE;
419   PetscBool              use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
420   PetscViewer            matl_dbg_viewer = NULL;
421   PetscErrorCode         ierr;
422 
423   PetscFunctionBegin;
424   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
425   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
426   if (Ain) {
427     ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
428     sub_schurs->A = Ain;
429   }
430 
431   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
432   sub_schurs->S = Sin;
433   if (sub_schurs->schur_explicit) {
434     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
435   }
436 
437   /* preliminary checks */
438   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");
439 
440   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
441 
442   /* debug (MATLAB) */
443   if (sub_schurs->debug) {
444     PetscMPIInt size,rank;
445     PetscInt    nr,*print_schurs_ranks,print_schurs = PETSC_FALSE;
446     PetscBool   flg;
447 
448     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);CHKERRQ(ierr);
449     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
450     nr   = size;
451     ierr = PetscMalloc1(nr,&print_schurs_ranks);CHKERRQ(ierr);
452     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
453     ierr = PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);CHKERRQ(ierr);
454     if (!flg) print_schurs = PETSC_TRUE;
455     else {
456       print_schurs = PETSC_FALSE;
457       for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
458     }
459     ierr = PetscOptionsEnd();CHKERRQ(ierr);
460     ierr = PetscFree(print_schurs_ranks);CHKERRQ(ierr);
461     if (print_schurs) {
462       char filename[256];
463 
464       ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
465       ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&matl_dbg_viewer);CHKERRQ(ierr);
466       ierr = PetscViewerPushFormat(matl_dbg_viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
467     }
468   }
469 
470 
471   /* restrict work on active processes */
472   if (sub_schurs->restrict_comm) {
473     PetscSubcomm subcomm;
474     PetscMPIInt  color,rank;
475 
476     color = 0;
477     if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
478     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
479     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
480     ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
481     ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
482     ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
483     ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
484     if (!sub_schurs->n_subs) {
485       ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
486       PetscFunctionReturn(0);
487     }
488   } else {
489     ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr);
490   }
491 
492   /* get Schur complement matrices */
493   if (!sub_schurs->schur_explicit) {
494     Mat       tA_IB,tA_BI,tA_BB;
495     PetscBool isseqsbaij;
496     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
497     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
498     if (isseqsbaij) {
499       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
500       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
501       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
502     } else {
503       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
504       A_BB = tA_BB;
505       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
506       A_IB = tA_IB;
507       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
508       A_BI = tA_BI;
509     }
510   } else {
511     A_II = NULL;
512     A_IB = NULL;
513     A_BI = NULL;
514     A_BB = NULL;
515   }
516   S_all = NULL;
517 
518   /* determine interior problems */
519   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
520   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
521     PetscBT                touched;
522     const PetscInt*        idx_B;
523     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
524 
525     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
526     /* get sizes */
527     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
528     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
529 
530     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
531     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
532     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
533 
534     /* all boundary dofs must be skipped when adding layers */
535     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
536     for (j=0;j<n_B;j++) {
537       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
538     }
539     ierr = PetscArraycpy(local_numbering,idx_B,n_B);CHKERRQ(ierr);
540     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
541 
542     /* add prescribed number of layers of dofs */
543     n_local_dofs = n_B;
544     n_prev_added = n_B;
545     for (layer=0;layer<nlayers;layer++) {
546       PetscInt n_added;
547       if (n_local_dofs == n_I+n_B) break;
548       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);
549       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
550       n_prev_added = n_added;
551       n_local_dofs += n_added;
552       if (!n_added) break;
553     }
554     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
555 
556     /* IS for I layer dofs in original numbering */
557     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);
558     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
559     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
560     /* IS for I layer dofs in I numbering */
561     if (!sub_schurs->schur_explicit) {
562       ISLocalToGlobalMapping ItoNmap;
563       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
564       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
565       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
566 
567       /* II block */
568       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
569     }
570   } else {
571     PetscInt n_I;
572 
573     /* IS for I dofs in original numbering */
574     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
575     is_I_layer = sub_schurs->is_I;
576 
577     /* IS for I dofs in I numbering (strided 1) */
578     if (!sub_schurs->schur_explicit) {
579       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
580       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
581 
582       /* II block is the same */
583       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
584       AE_II = A_II;
585     }
586   }
587 
588   /* Get info on subset sizes and sum of all subsets sizes */
589   max_subset_size = 0;
590   local_size = 0;
591   for (i=0;i<sub_schurs->n_subs;i++) {
592     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
593     max_subset_size = PetscMax(subset_size,max_subset_size);
594     local_size += subset_size;
595   }
596 
597   /* Work arrays for local indices */
598   extra = 0;
599   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
600   if (sub_schurs->schur_explicit && is_I_layer) {
601     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
602   }
603   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
604   if (extra) {
605     const PetscInt *idxs;
606     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
607     ierr = PetscArraycpy(all_local_idx_N,idxs,extra);CHKERRQ(ierr);
608     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
609   }
610   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
611   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
612   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
613 
614   /* Get local indices in local numbering */
615   local_size = 0;
616   local_stash_size = 0;
617   for (i=0;i<sub_schurs->n_subs;i++) {
618     PetscInt j;
619     const    PetscInt *idxs;
620 
621     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
622     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
623     /* start (smallest in global ordering) and multiplicity */
624     auxnum1[i] = idxs[0];
625     auxnum2[i] = subset_size*subset_size;
626     /* subset indices in local numbering */
627     ierr = PetscArraycpy(all_local_idx_N+local_size+extra,idxs,subset_size);CHKERRQ(ierr);
628     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
629     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
630     local_size += subset_size;
631     local_stash_size += subset_size*subset_size;
632   }
633 
634   /* allocate extra workspace needed only for GETRI or SYTRF */
635   use_potr = use_sytr = PETSC_FALSE;
636   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
637     use_potr = PETSC_TRUE;
638   } else if (sub_schurs->is_symmetric) {
639     use_sytr = PETSC_TRUE;
640   }
641   if (local_size && !use_potr) {
642     PetscScalar  lwork,dummyscalar = 0.;
643     PetscBLASInt dummyint = 0;
644 
645     B_lwork = -1;
646     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
647     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
648     if (use_sytr) {
649       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
650       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
651     } else {
652       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
653       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
654     }
655     ierr = PetscFPTrapPop();CHKERRQ(ierr);
656     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
657     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
658   } else {
659     Bwork = NULL;
660     pivots = NULL;
661   }
662 
663   /* prepare data for summing up properly schurs on subsets */
664   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
665   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
666   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
667   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
668   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
669   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
670   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
671   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
672   if (i != local_stash_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_stash_size);
673   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,local_stash_size,NULL,&lstash);CHKERRQ(ierr);
674   ierr = VecCreateMPI(comm_n,PETSC_DECIDE,global_size,&gstash);CHKERRQ(ierr);
675   ierr = VecScatterCreate(lstash,NULL,gstash,all_subsets_n,&sstash);CHKERRQ(ierr);
676   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
677 
678   /* subset indices in local boundary numbering */
679   if (!sub_schurs->is_Ej_all) {
680     PetscInt *all_local_idx_B;
681 
682     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
683     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
684     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
685     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
686   }
687 
688   if (change) {
689     ISLocalToGlobalMapping BtoS;
690     IS                     change_primal_B;
691     IS                     change_primal_all;
692 
693     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
694     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
695     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
696     for (i=0;i<sub_schurs->n_subs;i++) {
697       ISLocalToGlobalMapping NtoS;
698       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
699       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
700       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
701     }
702     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
703     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
704     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
705     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
706     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
707     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
708     for (i=0;i<sub_schurs->n_subs;i++) {
709       Mat change_sub;
710 
711       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
712       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
713       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
714       if (!sub_schurs->change_with_qr) {
715         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
716       } else {
717         Mat change_subt;
718         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
719         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
720         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
721       }
722       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
723       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
724       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
725       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
726     }
727     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
728   }
729 
730   /* Local matrix of all local Schur on subsets (transposed) */
731   if (!sub_schurs->S_Ej_all) {
732     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
733     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
734     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
735     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
736   }
737 
738   /* Compute Schur complements explicitly */
739   F = NULL;
740   if (!sub_schurs->schur_explicit) {
741     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
742        it is not efficient, unless the economic version of the scaling is used */
743     Mat         S_Ej_expl;
744     PetscScalar *work;
745     PetscInt    j,*dummy_idx;
746     PetscBool   Sdense;
747 
748     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
749     local_size = 0;
750     for (i=0;i<sub_schurs->n_subs;i++) {
751       IS  is_subset_B;
752       Mat AE_EE,AE_IE,AE_EI,S_Ej;
753 
754       /* subsets in original and boundary numbering */
755       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
756       /* EE block */
757       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
758       /* IE block */
759       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
760       /* EI block */
761       if (sub_schurs->is_symmetric) {
762         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
763       } else if (sub_schurs->is_hermitian) {
764         ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
765       } else {
766         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
767       }
768       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
769       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
770       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
771       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
772       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
773       if (AE_II == A_II) { /* we can reuse the same ksp */
774         KSP ksp;
775         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
776         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
777       } else { /* build new ksp object which inherits ksp and pc types from the original one */
778         KSP       origksp,schurksp;
779         PC        origpc,schurpc;
780         KSPType   ksp_type;
781         PetscInt  n_internal;
782         PetscBool ispcnone;
783 
784         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
785         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
786         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
787         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
788         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
789         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
790         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
791         if (!ispcnone) {
792           PCType pc_type;
793           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
794           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
795         } else {
796           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
797         }
798         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
799         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
800           MatSolverType solver = NULL;
801           ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr);
802           if (solver) {
803             ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr);
804           }
805         }
806         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
807       }
808       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
809       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
810       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
811       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
812       if (Sdense) {
813         for (j=0;j<subset_size;j++) {
814           dummy_idx[j]=local_size+j;
815         }
816         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
817       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
818       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
819       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
820       local_size += subset_size;
821     }
822     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
823     /* free */
824     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
825     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
826     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
827   } else {
828     Mat               A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
829     Vec               Dall = NULL;
830     IS                is_A_all,*is_p_r = NULL;
831     MatType           Stype;
832     PetscScalar       *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
833     const PetscScalar *rS_data;
834     PetscInt          n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
835     PetscBool         economic,solver_S,S_lower_triangular = PETSC_FALSE;
836     PetscBool         schur_has_vertices,factor_workaround;
837     PetscBool         use_cholesky;
838 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
839     PetscBool         oldpin;
840 #endif
841 
842     /* get sizes */
843     n_I = 0;
844     if (is_I_layer) {
845       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
846     }
847     economic = PETSC_FALSE;
848     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
849     if (cum != n_I) economic = PETSC_TRUE;
850     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
851     size_active_schur = local_size;
852 
853     /* import scaling vector (wrong formulation if we have 3D edges) */
854     if (scaling && compute_Stilda) {
855       const PetscScalar *array;
856       PetscScalar       *array2;
857       const PetscInt    *idxs;
858       PetscInt          i;
859 
860       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
861       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
862       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
863       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
864       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
865       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
866       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
867       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
868       deluxe = PETSC_FALSE;
869     }
870 
871     /* size active schurs does not count any dirichlet or vertex dof on the interface */
872     factor_workaround = PETSC_FALSE;
873     schur_has_vertices = PETSC_FALSE;
874     cum = n_I+size_active_schur;
875     if (sub_schurs->is_dir) {
876       const PetscInt* idxs;
877       PetscInt        n_dir;
878 
879       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
880       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
881       ierr = PetscArraycpy(all_local_idx_N+cum,idxs,n_dir);CHKERRQ(ierr);
882       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
883       cum += n_dir;
884       factor_workaround = PETSC_TRUE;
885     }
886     /* include the primal vertices in the Schur complement */
887     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
888       PetscInt n_v;
889 
890       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
891       if (n_v) {
892         const PetscInt* idxs;
893 
894         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
895         ierr = PetscArraycpy(all_local_idx_N+cum,idxs,n_v);CHKERRQ(ierr);
896         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
897         cum += n_v;
898         factor_workaround = PETSC_TRUE;
899         schur_has_vertices = PETSC_TRUE;
900       }
901     }
902     size_schur = cum - n_I;
903     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
904 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
905     oldpin = sub_schurs->A->pinnedtocpu;
906     ierr = MatPinToCPU(sub_schurs->A,PETSC_TRUE);CHKERRQ(ierr);
907 #endif
908     if (cum == n) {
909       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
910       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
911     } else {
912       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
913     }
914 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
915     ierr = MatPinToCPU(sub_schurs->A,oldpin);CHKERRQ(ierr);
916 #endif
917     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
918     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
919 
920     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
921        this is a workaround */
922     if (benign_n) {
923       Vec                    v,benign_AIIm1_ones;
924       ISLocalToGlobalMapping N_to_reor;
925       IS                     is_p0,is_p0_p;
926       PetscScalar            *cs_AIB,*AIIm1_data;
927       PetscInt               sizeA;
928 
929       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
930       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
931       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
932       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
933       ierr = MatCreateVecs(A,&v,&benign_AIIm1_ones);CHKERRQ(ierr);
934       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
935       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
936       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
937       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
938       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
939       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
940       /* compute colsum of A_IB restricted to pressures */
941       for (i=0;i<benign_n;i++) {
942         const PetscScalar *array;
943         const PetscInt    *idxs;
944         PetscInt          j,nz;
945 
946         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
947         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
948         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
949         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
950         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
951         ierr = VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);CHKERRQ(ierr);
952         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
953         ierr = VecResetArray(benign_AIIm1_ones);CHKERRQ(ierr);
954         ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
955         for (j=0;j<size_schur;j++) {
956 #if defined(PETSC_USE_COMPLEX)
957           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
958 #else
959           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
960 #endif
961         }
962         ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
963       }
964       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
965       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
966       ierr = VecDestroy(&v);CHKERRQ(ierr);
967       ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
968       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
969       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
970       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
971       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
972       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
973       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
974     }
975     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr);
976     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
977     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
978 
979     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
980     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
981 
982     /* when using the benign subspace trick, the local Schur complements are SPD */
983     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
984 
985     if (n_I) {
986       IS        is_schur;
987       char      stype[64];
988       PetscBool gpu;
989 
990       if (use_cholesky) {
991         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
992       } else {
993         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
994       }
995       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
996 
997       /* subsets ordered last */
998       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
999       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
1000       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
1001 
1002       /* factorization step */
1003       if (use_cholesky) {
1004         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
1005 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1006         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
1007 #endif
1008         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1009         S_lower_triangular = PETSC_TRUE;
1010       } else {
1011         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
1012 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1013         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
1014 #endif
1015         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1016       }
1017       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
1018 
1019       if (matl_dbg_viewer) {
1020         Mat S;
1021         IS  is;
1022 
1023         ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
1024         ierr = MatView(A,matl_dbg_viewer);CHKERRQ(ierr);
1025         ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
1026         ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr);
1027         ierr = MatView(S,matl_dbg_viewer);CHKERRQ(ierr);
1028         ierr = MatDestroy(&S);CHKERRQ(ierr);
1029         ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr);
1030         ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr);
1031         ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr);
1032         ierr = ISDestroy(&is);CHKERRQ(ierr);
1033         ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr);
1034         ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr);
1035         ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr);
1036         ierr = ISDestroy(&is);CHKERRQ(ierr);
1037         ierr = PetscObjectSetName((PetscObject)is_A_all,"IA");CHKERRQ(ierr);
1038         ierr = ISView(is_A_all,matl_dbg_viewer);CHKERRQ(ierr);
1039       }
1040 
1041       /* get explicit Schur Complement computed during numeric factorization */
1042       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1043       ierr = PetscStrncpy(stype,MATSEQDENSE,sizeof(stype));CHKERRQ(ierr);
1044       ierr = PetscObjectTypeCompareAny((PetscObject)A,&gpu,MATSEQAIJVIENNACL,MATSEQAIJCUSPARSE,"");CHKERRQ(ierr);
1045       if (gpu) {
1046         ierr = PetscStrncpy(stype,MATSEQDENSECUDA,sizeof(stype));CHKERRQ(ierr);
1047       }
1048       ierr = PetscOptionsGetString(NULL,sub_schurs->prefix,"-sub_schurs_schur_mat_type",stype,sizeof(stype),NULL);CHKERRQ(ierr);
1049       ierr = MatConvert(S_all,stype,MAT_INPLACE_MATRIX,&S_all);CHKERRQ(ierr);
1050       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
1051       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
1052       ierr = MatGetType(S_all,&Stype);CHKERRQ(ierr);
1053 
1054       /* we can reuse the solvers if we are not using the economic version */
1055       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1056       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1057       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
1058         reuse_solvers = factor_workaround = PETSC_FALSE;
1059 
1060       solver_S = PETSC_TRUE;
1061 
1062       /* update the Schur complement with the change of basis on the pressures */
1063       if (benign_n) {
1064         const PetscScalar *cs_AIB;
1065         PetscScalar       *S_data,*AIIm1_data;
1066         Mat               S2 = NULL,S3 = NULL; /* dbg */
1067         PetscScalar       *S2_data,*S3_data; /* dbg */
1068         Vec               v,benign_AIIm1_ones;
1069         PetscInt          sizeA;
1070 
1071         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1072         ierr = MatCreateVecs(A,&v,&benign_AIIm1_ones);CHKERRQ(ierr);
1073         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1074 #if defined(PETSC_HAVE_MUMPS)
1075         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1076 #endif
1077 #if defined(PETSC_HAVE_MKL_PARDISO)
1078         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1079 #endif
1080         ierr = MatDenseGetArrayRead(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1081         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1082         if (matl_dbg_viewer) {
1083           ierr = MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S2);CHKERRQ(ierr);
1084           ierr = MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S3);CHKERRQ(ierr);
1085           ierr = MatDenseGetArray(S2,&S2_data);CHKERRQ(ierr);
1086           ierr = MatDenseGetArray(S3,&S3_data);CHKERRQ(ierr);
1087         }
1088         for (i=0;i<benign_n;i++) {
1089           PetscScalar    *array,sum = 0.,one = 1.,*sums;
1090           const PetscInt *idxs;
1091           PetscInt       k,j,nz;
1092           PetscBLASInt   B_k,B_n;
1093 
1094           ierr = PetscCalloc1(benign_n,&sums);CHKERRQ(ierr);
1095           ierr = VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);CHKERRQ(ierr);
1096           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
1097           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1098           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1099           ierr = VecResetArray(benign_AIIm1_ones);CHKERRQ(ierr);
1100           /* p0 dofs (eliminated) are excluded from the sums */
1101           for (k=0;k<benign_n;k++) {
1102             ierr = ISGetLocalSize(is_p_r[k],&nz);CHKERRQ(ierr);
1103             ierr = ISGetIndices(is_p_r[k],&idxs);CHKERRQ(ierr);
1104             for (j=0;j<nz-1;j++) sums[k] -= AIIm1_data[idxs[j]+sizeA*i];
1105             ierr = ISRestoreIndices(is_p_r[k],&idxs);CHKERRQ(ierr);
1106           }
1107           ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
1108           if (matl_dbg_viewer) {
1109             Vec  vv;
1110             char name[16];
1111 
1112             ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array+n_I,&vv);CHKERRQ(ierr);
1113             ierr = PetscSNPrintf(name,sizeof(name),"Pvs%D",i);CHKERRQ(ierr);
1114             ierr = PetscObjectSetName((PetscObject)vv,name);CHKERRQ(ierr);
1115             ierr = VecView(vv,matl_dbg_viewer);CHKERRQ(ierr);
1116           }
1117           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1118           /* cs_AIB already scaled by 1./nz */
1119           B_k = 1;
1120           for (k=0;k<benign_n;k++) {
1121             sum  = sums[k];
1122             ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
1123 
1124             if (PetscAbsScalar(sum) == 0.0) continue;
1125             if (k == i) {
1126               PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1127               if (matl_dbg_viewer) {
1128                 PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n));
1129               }
1130             } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */
1131               sum /= 2.0;
1132               PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1133               if (matl_dbg_viewer) {
1134                 PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n));
1135               }
1136             }
1137           }
1138           sum  = 1.;
1139           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));
1140           if (matl_dbg_viewer) {
1141             PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S2_data,&B_n));
1142           }
1143           ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
1144           /* set p0 entry of AIIm1_ones to zero */
1145           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1146           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1147           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1148           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1149           ierr = PetscFree(sums);CHKERRQ(ierr);
1150         }
1151         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1152         if (matl_dbg_viewer) {
1153           ierr = MatDenseRestoreArray(S2,&S2_data);CHKERRQ(ierr);
1154           ierr = MatDenseRestoreArray(S3,&S3_data);CHKERRQ(ierr);
1155         }
1156         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1157           PetscInt k,j;
1158           for (k=0;k<size_schur;k++) {
1159             for (j=k;j<size_schur;j++) {
1160               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1161             }
1162           }
1163         }
1164 
1165         /* restore defaults */
1166 #if defined(PETSC_HAVE_MUMPS)
1167         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
1168 #endif
1169 #if defined(PETSC_HAVE_MKL_PARDISO)
1170         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
1171 #endif
1172         ierr = MatDenseRestoreArrayRead(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1173         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1174         ierr = VecDestroy(&v);CHKERRQ(ierr);
1175         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1176         if (matl_dbg_viewer) {
1177           Mat S;
1178 
1179           ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1180           ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
1181           ierr = PetscObjectSetName((PetscObject)S,"Sb");CHKERRQ(ierr);
1182           ierr = MatView(S,matl_dbg_viewer);CHKERRQ(ierr);
1183           ierr = MatDestroy(&S);CHKERRQ(ierr);
1184           ierr = PetscObjectSetName((PetscObject)S2,"S2P");CHKERRQ(ierr);
1185           ierr = MatView(S2,matl_dbg_viewer);CHKERRQ(ierr);
1186           ierr = PetscObjectSetName((PetscObject)S3,"S3P");CHKERRQ(ierr);
1187           ierr = MatView(S3,matl_dbg_viewer);CHKERRQ(ierr);
1188           ierr = PetscObjectSetName((PetscObject)cs_AIB_mat,"cs");CHKERRQ(ierr);
1189           ierr = MatView(cs_AIB_mat,matl_dbg_viewer);CHKERRQ(ierr);
1190           ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1191         }
1192         ierr = MatDestroy(&S2);CHKERRQ(ierr);
1193         ierr = MatDestroy(&S3);CHKERRQ(ierr);
1194       }
1195       if (!reuse_solvers) {
1196         for (i=0;i<benign_n;i++) {
1197           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1198         }
1199         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
1200         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
1201         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1202       }
1203     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1204       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1205       ierr = MatGetType(S_all,&Stype);CHKERRQ(ierr);
1206       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1207       factor_workaround = PETSC_FALSE;
1208       solver_S = PETSC_FALSE;
1209     }
1210 
1211     if (reuse_solvers) {
1212       Mat                A_II,Afake;
1213       Vec                vec1_B;
1214       PCBDDCReuseSolvers msolv_ctx;
1215       PetscInt           n_R;
1216 
1217       if (sub_schurs->reuse_solver) {
1218         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1219       } else {
1220         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1221       }
1222       msolv_ctx = sub_schurs->reuse_solver;
1223       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1224       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1225       msolv_ctx->F = F;
1226       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1227       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1228       {
1229         PetscScalar *array;
1230         PetscInt    n;
1231 
1232         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1233         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1234         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1235         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1236       }
1237       msolv_ctx->has_vertices = schur_has_vertices;
1238 
1239       /* interior solver */
1240       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1241       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1242       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1243       ierr = PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");CHKERRQ(ierr);
1244       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1245       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1246       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1247       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1248 
1249       /* correction solver */
1250       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1251       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1252       ierr = PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");CHKERRQ(ierr);
1253       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1254       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1255       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1256       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
1257 
1258       /* scatter and vecs for Schur complement solver */
1259       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
1260       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1261       if (!schur_has_vertices) {
1262         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1263         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1264         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
1265         msolv_ctx->is_R = is_A_all;
1266       } else {
1267         IS              is_B_all;
1268         const PetscInt* idxs;
1269         PetscInt        dual,n_v,n;
1270 
1271         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1272         dual = size_schur - n_v;
1273         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1274         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1275         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1276         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1277         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1278         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1279         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1280         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1281         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1282         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1283       }
1284       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
1285       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
1286       ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1287       ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1288       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
1289       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1290       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1291 
1292       /* communicate benign info to solver context */
1293       if (benign_n) {
1294         PetscScalar *array;
1295 
1296         msolv_ctx->benign_n = benign_n;
1297         msolv_ctx->benign_zerodiag_subs = is_p_r;
1298         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
1299         msolv_ctx->benign_csAIB = cs_AIB_mat;
1300         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
1301         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1302         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1303         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1304         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1305       }
1306     } else {
1307       if (sub_schurs->reuse_solver) {
1308         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1309       }
1310       ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1311     }
1312     ierr = MatDestroy(&A);CHKERRQ(ierr);
1313     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
1314 
1315     /* Work arrays */
1316     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
1317 
1318     /* matrices for deluxe scaling and adaptive selection */
1319     if (compute_Stilda) {
1320       if (!sub_schurs->sum_S_Ej_tilda_all) {
1321         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1322       }
1323       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1324         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1325       }
1326     }
1327 
1328     /* S_Ej_all */
1329     cum = cum2 = 0;
1330     ierr = MatDenseGetArrayRead(S_all,&rS_data);CHKERRQ(ierr);
1331     for (i=0;i<sub_schurs->n_subs;i++) {
1332       PetscInt j;
1333 
1334       /* get S_E */
1335       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1336       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1337         PetscInt k;
1338         for (k=0;k<subset_size;k++) {
1339           for (j=k;j<subset_size;j++) {
1340             work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1341             work[j*subset_size+k] = PetscConj(rS_data[cum2+k*size_schur+j]);
1342           }
1343         }
1344       } else { /* just copy to workspace */
1345         PetscInt k;
1346         for (k=0;k<subset_size;k++) {
1347           for (j=0;j<subset_size;j++) {
1348             work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1349           }
1350         }
1351       }
1352       /* insert S_E values */
1353       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1354       if (sub_schurs->change) {
1355         Mat change_sub,SEj,T;
1356 
1357         /* change basis */
1358         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1359         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1360         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1361           Mat T2;
1362           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1363           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1364           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1365           ierr = MatDestroy(&T2);CHKERRQ(ierr);
1366         } else {
1367           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1368         }
1369         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1370         ierr = MatDestroy(&T);CHKERRQ(ierr);
1371         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
1372         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1373       }
1374       if (deluxe) {
1375         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1376         /* if adaptivity is requested, invert S_E blocks */
1377         if (compute_Stilda) {
1378           Mat               M;
1379           const PetscScalar *vals;
1380           PetscBool         isdense,isdensecuda;
1381 
1382           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&M);CHKERRQ(ierr);
1383           ierr = MatSetOption(M,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
1384           ierr = MatSetOption(M,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
1385           if (!PetscBTLookup(sub_schurs->is_edge,i)) {
1386             ierr = MatSetType(M,Stype);CHKERRQ(ierr);
1387           }
1388           ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1389           ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQDENSECUDA,&isdensecuda);CHKERRQ(ierr);
1390           if (use_cholesky) {
1391             ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
1392           } else {
1393             ierr = MatLUFactor(M,NULL,NULL,NULL);CHKERRQ(ierr);
1394           }
1395           if (isdense) {
1396             ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
1397 #if defined(PETSC_HAVE_CUDA)
1398           } else if (isdensecuda) {
1399             ierr = MatSeqDenseCUDAInvertFactors_Private(M);CHKERRQ(ierr);
1400 #endif
1401           } else SETERRQ1(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"Not implemented for type %s",Stype);
1402           ierr = MatDenseGetArrayRead(M,&vals);CHKERRQ(ierr);
1403           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
1404           ierr = MatDenseRestoreArrayRead(M,&vals);CHKERRQ(ierr);
1405           ierr = MatDestroy(&M);CHKERRQ(ierr);
1406         }
1407       } else if (compute_Stilda) { /* not using deluxe */
1408         Mat         SEj;
1409         Vec         D;
1410         PetscScalar *array;
1411 
1412         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1413         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
1414         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
1415         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1416         ierr = VecShift(D,-1.);CHKERRQ(ierr);
1417         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
1418         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1419         ierr = VecDestroy(&D);CHKERRQ(ierr);
1420         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1421       }
1422       cum += subset_size;
1423       cum2 += subset_size*(size_schur + 1);
1424     }
1425     ierr = MatDenseRestoreArrayRead(S_all,&rS_data);CHKERRQ(ierr);
1426 
1427     if (solver_S) {
1428       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1429     }
1430 
1431     /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory
1432        however, preliminary tests indicate using GPUs is still faster in the solve phase */
1433 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1434     if (reuse_solvers) {
1435       Mat                  St;
1436       MatFactorSchurStatus st;
1437       PetscBool            flg = PETSC_FALSE;
1438 
1439       ierr = PetscOptionsGetBool(NULL,sub_schurs->prefix,"-sub_schurs_schur_pin_to_cpu",&flg,NULL);CHKERRQ(ierr);
1440       ierr = MatFactorGetSchurComplement(F,&St,&st);CHKERRQ(ierr);
1441       ierr = MatPinToCPU(St,flg);CHKERRQ(ierr);
1442       ierr = MatFactorRestoreSchurComplement(F,&St,st);CHKERRQ(ierr);
1443     }
1444 #endif
1445 
1446     schur_factor = NULL;
1447     if (compute_Stilda && size_active_schur) {
1448 
1449       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1450         PetscInt j;
1451 
1452         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1453         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1454       } else {
1455         Mat S_all_inv=NULL;
1456 
1457         if (solver_S) {
1458           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1459              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 */
1460           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1461             PetscScalar *data;
1462             PetscInt     nd = 0;
1463 
1464             if (!use_potr) {
1465               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1466             }
1467             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1468             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1469             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1470               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1471             }
1472 
1473             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1474             if (schur_has_vertices) {
1475               Mat          M;
1476               PetscScalar *tdata;
1477               PetscInt     nv = 0, news;
1478 
1479               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1480               news = size_active_schur + nv;
1481               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1482               for (i=0;i<size_active_schur;i++) {
1483                 ierr = PetscArraycpy(tdata+i*(news+1),data+i*(size_schur+1),size_active_schur-i);CHKERRQ(ierr);
1484                 ierr = PetscArraycpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv);CHKERRQ(ierr);
1485               }
1486               for (i=0;i<nv;i++) {
1487                 PetscInt k = i+size_active_schur;
1488                 ierr = PetscArraycpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),nv-i);CHKERRQ(ierr);
1489               }
1490 
1491               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
1492               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1493               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
1494               /* save the factors */
1495               cum = 0;
1496               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
1497               for (i=0;i<size_active_schur;i++) {
1498                 ierr = PetscArraycpy(schur_factor+cum,tdata+i*(news+1),size_active_schur-i);CHKERRQ(ierr);
1499                 cum += size_active_schur - i;
1500               }
1501               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1502               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
1503               /* move back just the active dofs to the Schur complement */
1504               for (i=0;i<size_active_schur;i++) {
1505                 ierr = PetscArraycpy(data+i*size_schur,tdata+i*news,size_active_schur);CHKERRQ(ierr);
1506               }
1507               ierr = PetscFree(tdata);CHKERRQ(ierr);
1508               ierr = MatDestroy(&M);CHKERRQ(ierr);
1509             } else { /* we can factorize and invert just the activedofs */
1510               Mat         M;
1511               PetscScalar *aux;
1512 
1513               ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr);
1514               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1515               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
1516               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
1517               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1518               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
1519               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
1520               ierr = MatDestroy(&M);CHKERRQ(ierr);
1521               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr);
1522               ierr = MatZeroEntries(M);CHKERRQ(ierr);
1523               ierr = MatDestroy(&M);CHKERRQ(ierr);
1524               ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr);
1525               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
1526               ierr = MatZeroEntries(M);CHKERRQ(ierr);
1527               ierr = MatDestroy(&M);CHKERRQ(ierr);
1528               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1529               ierr = PetscFree(aux);CHKERRQ(ierr);
1530             }
1531             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1532           } else { /* use MatFactor calls to invert S */
1533             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
1534             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1535           }
1536         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1537           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1538           S_all_inv = S_all;
1539           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1540           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1541           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1542           if (use_potr) {
1543             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1544             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1545             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1546             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1547           } else if (use_sytr) {
1548             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1549             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1550             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1551             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1552           } else {
1553             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1554             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1555             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1556             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1557           }
1558           ierr = PetscLogFlops(1.0*size_schur*size_schur*size_schur);CHKERRQ(ierr);
1559           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1560           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1561         }
1562         /* S_Ej_tilda_all */
1563         cum = cum2 = 0;
1564         ierr = MatDenseGetArrayRead(S_all_inv,&rS_data);CHKERRQ(ierr);
1565         for (i=0;i<sub_schurs->n_subs;i++) {
1566           PetscInt j;
1567 
1568           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1569           /* get (St^-1)_E */
1570           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1571              will be properly accessed later during adaptive selection */
1572           if (S_lower_triangular) {
1573             PetscInt k;
1574             if (sub_schurs->change) {
1575               for (k=0;k<subset_size;k++) {
1576                 for (j=k;j<subset_size;j++) {
1577                   work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1578                   work[j*subset_size+k] = work[k*subset_size+j];
1579                 }
1580               }
1581             } else {
1582               for (k=0;k<subset_size;k++) {
1583                 for (j=k;j<subset_size;j++) {
1584                   work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1585                 }
1586               }
1587             }
1588           } else {
1589             PetscInt k;
1590             for (k=0;k<subset_size;k++) {
1591               for (j=0;j<subset_size;j++) {
1592                 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1593               }
1594             }
1595           }
1596           if (sub_schurs->change) {
1597             Mat change_sub,SEj,T;
1598 
1599             /* change basis */
1600             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1601             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1602             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1603               Mat T2;
1604               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1605               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1606               ierr = MatDestroy(&T2);CHKERRQ(ierr);
1607               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1608             } else {
1609               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1610             }
1611             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1612             ierr = MatDestroy(&T);CHKERRQ(ierr);
1613             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1614             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
1615             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1616           }
1617           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1618           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1619           cum += subset_size;
1620           cum2 += subset_size*(size_schur + 1);
1621         }
1622         ierr = MatDenseRestoreArrayRead(S_all_inv,&rS_data);CHKERRQ(ierr);
1623         if (solver_S) {
1624           if (schur_has_vertices) {
1625             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
1626           } else {
1627             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
1628           }
1629         }
1630         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1631       }
1632 
1633       /* move back factors if needed */
1634       if (schur_has_vertices) {
1635         Mat      S_tmp;
1636         PetscInt nd = 0;
1637 
1638         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1639         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1640         if (use_potr) {
1641           PetscScalar *data;
1642 
1643           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1644           ierr = PetscArrayzero(data,size_schur*size_schur);CHKERRQ(ierr);
1645 
1646           if (S_lower_triangular) {
1647             cum = 0;
1648             for (i=0;i<size_active_schur;i++) {
1649               ierr = PetscArraycpy(data+i*(size_schur+1),schur_factor+cum,size_active_schur-i);CHKERRQ(ierr);
1650               cum += size_active_schur-i;
1651 	    }
1652           } else {
1653             ierr = PetscArraycpy(data,schur_factor,size_schur*size_schur);CHKERRQ(ierr);
1654           }
1655           if (sub_schurs->is_dir) {
1656             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1657 	    for (i=0;i<nd;i++) {
1658 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1659 	    }
1660           }
1661           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1662              set the diagonal entry of the Schur factor to a very large value */
1663           for (i=size_active_schur+nd;i<size_schur;i++) {
1664             data[i*(size_schur+1)] = infty;
1665           }
1666           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1667         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1668         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
1669       }
1670     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1671       PetscScalar *data;
1672       PetscInt    nd = 0;
1673 
1674       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1675         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1676       }
1677       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1678       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1679       for (i=0;i<size_active_schur;i++) {
1680         ierr = PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);CHKERRQ(ierr);
1681       }
1682       for (i=size_active_schur+nd;i<size_schur;i++) {
1683         ierr = PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);CHKERRQ(ierr);
1684         data[i*(size_schur+1)] = infty;
1685       }
1686       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1687       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1688     }
1689     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1690     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1691     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
1692   }
1693   ierr = PetscFree(nnz);CHKERRQ(ierr);
1694   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1695   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1696   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1697   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1698   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1699   ierr = MatDestroy(&F);CHKERRQ(ierr);
1700 
1701   /* speed up matrix assembly */
1702   ierr = PetscMalloc1(sub_schurs->n_subs,&nnz);CHKERRQ(ierr);
1703   for (i=0;i<sub_schurs->n_subs;i++) {
1704     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);CHKERRQ(ierr);
1705   }
1706   ierr = ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);CHKERRQ(ierr);
1707   ierr = MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1708   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1709   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1710   if (compute_Stilda) {
1711     ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1712     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1713     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1714     if (deluxe) {
1715       ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1716       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1717       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1718     }
1719   }
1720   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1721 
1722   /* Get local part of (\sum_j S_Ej) */
1723   if (!sub_schurs->sum_S_Ej_all) {
1724     ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1725   }
1726   ierr = VecSet(gstash,0.0);CHKERRQ(ierr);
1727   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&stasharray);CHKERRQ(ierr);
1728   ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr);
1729   ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1730   ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1731   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&stasharray);CHKERRQ(ierr);
1732   ierr = VecResetArray(lstash);CHKERRQ(ierr);
1733   ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&stasharray);CHKERRQ(ierr);
1734   ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr);
1735   ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1736   ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1737   ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&stasharray);CHKERRQ(ierr);
1738   ierr = VecResetArray(lstash);CHKERRQ(ierr);
1739 
1740   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1741   if (compute_Stilda) {
1742     ierr = VecSet(gstash,0.0);CHKERRQ(ierr);
1743     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);CHKERRQ(ierr);
1744     ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr);
1745     ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1746     ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1747     ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1748     ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1749     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);CHKERRQ(ierr);
1750     ierr = VecResetArray(lstash);CHKERRQ(ierr);
1751     if (deluxe) {
1752       ierr = VecSet(gstash,0.0);CHKERRQ(ierr);
1753       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);CHKERRQ(ierr);
1754       ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr);
1755       ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1756       ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1757       ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1758       ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1759       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);CHKERRQ(ierr);
1760       ierr = VecResetArray(lstash);CHKERRQ(ierr);
1761     } else {
1762       PetscScalar *array;
1763       PetscInt    cum;
1764 
1765       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1766       cum = 0;
1767       for (i=0;i<sub_schurs->n_subs;i++) {
1768         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1769         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1770         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1771         if (use_potr) {
1772           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1773           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1774           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1775           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1776         } else if (use_sytr) {
1777           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1778           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1779           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1780           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1781         } else {
1782           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1783           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1784           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1785           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1786         }
1787         ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr);
1788         ierr = PetscFPTrapPop();CHKERRQ(ierr);
1789         cum += subset_size*subset_size;
1790       }
1791       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1792       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1793       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1794       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1795     }
1796   }
1797   ierr = VecDestroy(&lstash);CHKERRQ(ierr);
1798   ierr = VecDestroy(&gstash);CHKERRQ(ierr);
1799   ierr = VecScatterDestroy(&sstash);CHKERRQ(ierr);
1800 
1801   if (matl_dbg_viewer) {
1802     PetscInt cum;
1803 
1804     if (sub_schurs->S_Ej_all) {
1805       ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr);
1806       ierr = MatView(sub_schurs->S_Ej_all,matl_dbg_viewer);CHKERRQ(ierr);
1807     }
1808     if (sub_schurs->sum_S_Ej_all) {
1809       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr);
1810       ierr = MatView(sub_schurs->sum_S_Ej_all,matl_dbg_viewer);CHKERRQ(ierr);
1811     }
1812     if (sub_schurs->sum_S_Ej_inv_all) {
1813       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr);
1814       ierr = MatView(sub_schurs->sum_S_Ej_inv_all,matl_dbg_viewer);CHKERRQ(ierr);
1815     }
1816     if (sub_schurs->sum_S_Ej_tilda_all) {
1817       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr);
1818       ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,matl_dbg_viewer);CHKERRQ(ierr);
1819     }
1820     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1821       IS   is;
1822       char name[16];
1823 
1824       ierr = PetscSNPrintf(name,sizeof(name),"IE%D",i);CHKERRQ(ierr);
1825       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1826       ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr);
1827       ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr);
1828       ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr);
1829       ierr = ISDestroy(&is);CHKERRQ(ierr);
1830       cum += subset_size;
1831     }
1832   }
1833 
1834   /* free workspace */
1835   ierr = PetscViewerDestroy(&matl_dbg_viewer);CHKERRQ(ierr);
1836   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1837   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1838   PetscFunctionReturn(0);
1839 }
1840 
1841 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1842 {
1843   IS              *faces,*edges,*all_cc,vertices;
1844   PetscInt        i,n_faces,n_edges,n_all_cc;
1845   PetscBool       is_sorted,ispetsc;
1846   PetscErrorCode  ierr;
1847 
1848   PetscFunctionBegin;
1849   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1850   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1851   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1852   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1853 
1854   /* reset any previous data */
1855   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1856 
1857   /* get index sets for faces and edges (already sorted by global ordering) */
1858   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1859   n_all_cc = n_faces+n_edges;
1860   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1861   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1862   for (i=0;i<n_faces;i++) {
1863     if (copycc) {
1864       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
1865     } else {
1866       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1867       all_cc[i] = faces[i];
1868     }
1869   }
1870   for (i=0;i<n_edges;i++) {
1871     if (copycc) {
1872       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
1873     } else {
1874       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1875       all_cc[n_faces+i] = edges[i];
1876     }
1877     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1878   }
1879   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1880   sub_schurs->is_vertices = vertices;
1881   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1882   sub_schurs->is_dir = NULL;
1883   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1884 
1885   /* Determine if MatFactor can be used */
1886   ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr);
1887 #if defined(PETSC_HAVE_MUMPS)
1888   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);CHKERRQ(ierr);
1889 #elif defined(PETSC_HAVE_MKL_PARDISO)
1890   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);CHKERRQ(ierr);
1891 #else
1892   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);CHKERRQ(ierr);
1893 #endif
1894 #if defined(PETSC_USE_COMPLEX)
1895   sub_schurs->is_hermitian  = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1896 #else
1897   sub_schurs->is_hermitian  = PETSC_TRUE;
1898 #endif
1899   sub_schurs->is_posdef     = PETSC_TRUE;
1900   sub_schurs->is_symmetric  = PETSC_TRUE;
1901   sub_schurs->debug         = PETSC_FALSE;
1902   sub_schurs->restrict_comm = PETSC_FALSE;
1903   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
1904   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);
1905   ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr);
1906   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
1907   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
1908   ierr = PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);CHKERRQ(ierr);
1909   ierr = PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);CHKERRQ(ierr);
1910   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1911   ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);CHKERRQ(ierr);
1912   sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1913 
1914   /* for reals, symmetric and hermitian are synonims */
1915 #if !defined(PETSC_USE_COMPLEX)
1916   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1917   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1918 #endif
1919 
1920   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1921   sub_schurs->is_I = is_I;
1922   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1923   sub_schurs->is_B = is_B;
1924   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1925   sub_schurs->l2gmap = graph->l2gmap;
1926   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1927   sub_schurs->BtoNmap = BtoNmap;
1928   sub_schurs->n_subs = n_all_cc;
1929   sub_schurs->is_subs = all_cc;
1930   sub_schurs->S_Ej_all = NULL;
1931   sub_schurs->sum_S_Ej_all = NULL;
1932   sub_schurs->sum_S_Ej_inv_all = NULL;
1933   sub_schurs->sum_S_Ej_tilda_all = NULL;
1934   sub_schurs->is_Ej_all = NULL;
1935   PetscFunctionReturn(0);
1936 }
1937 
1938 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1939 {
1940   PCBDDCSubSchurs schurs_ctx;
1941   PetscErrorCode  ierr;
1942 
1943   PetscFunctionBegin;
1944   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1945   schurs_ctx->n_subs = 0;
1946   *sub_schurs = schurs_ctx;
1947   PetscFunctionReturn(0);
1948 }
1949 
1950 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1951 {
1952   PetscInt       i;
1953   PetscErrorCode ierr;
1954 
1955   PetscFunctionBegin;
1956   if (!sub_schurs) PetscFunctionReturn(0);
1957   ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr);
1958   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1959   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1960   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1961   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1962   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1963   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1964   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1965   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1966   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1967   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1968   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1969   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1970   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1971   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1972   for (i=0;i<sub_schurs->n_subs;i++) {
1973     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1974   }
1975   if (sub_schurs->n_subs) {
1976     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1977   }
1978   if (sub_schurs->reuse_solver) {
1979     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1980   }
1981   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1982   if (sub_schurs->change) {
1983     for (i=0;i<sub_schurs->n_subs;i++) {
1984       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1985       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
1986     }
1987   }
1988   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1989   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
1990   sub_schurs->n_subs = 0;
1991   PetscFunctionReturn(0);
1992 }
1993 
1994 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1995 {
1996   PetscErrorCode ierr;
1997 
1998   PetscFunctionBegin;
1999   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
2000   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
2001   PetscFunctionReturn(0);
2002 }
2003 
2004 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
2005 {
2006   PetscInt       i,j,n;
2007   PetscErrorCode ierr;
2008 
2009   PetscFunctionBegin;
2010   n = 0;
2011   for (i=-n_prev;i<0;i++) {
2012     PetscInt start_dof = queue_tip[i];
2013     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
2014       PetscInt dof = adjncy[j];
2015       if (!PetscBTLookup(touched,dof)) {
2016         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
2017         queue_tip[n] = dof;
2018         n++;
2019       }
2020     }
2021   }
2022   *n_added = n;
2023   PetscFunctionReturn(0);
2024 }
2025