[petsc-users] Question about SuperLU

Xiaoye S. Li xsli at lbl.gov
Fri Jun 10 19:44:05 CDT 2022


On Fri, Jun 10, 2022 at 7:43 AM Tang, Qi <tangqi at msu.edu> wrote:

> We use superlu_dist.
>
> We have a 2 x 2 block where directly calling suplerlu_dist fails, but a pc
> based on a fieldsplit Schur complement + superlu_dist on the assembled
> Schur complement matrix converges. (All the converge criteria are default
> at this level)
>
> I am having a hard time to understand what is going on. The B,V block is
> of size 240K, so it is also hard to analyze. And the mat is not something
> we explicitly formed. It is formed by finite difference coloring Jacobian +
> a few levels of Schur complement.
>
>   / A 0 \
>   \ 0  I /
>
> Matt, I do not see this can explain why the second pc with superlu on S =
> A would succeed, if A is not full rank.
>
> I believe I found somewhere it says petsc’s pclu (or maybe superlu_dist)
> did reordering and it may introduce 0 pivoting. We are asking because it
> seems there is something we do not understand from pclu/superlu level.
>

If the matrix is non-singular, and you use the RowPerm default option:
LargeDiag_MC64, then you won't have zero pivot (in a structural sense),
unless numerical cancellation causes the diagonal element underflow, then
flush to zero.

You can try to set ReplaceTinyPivot:
-mat_superlu_dist_replacetinypivot

replace tiny pivots
See my reply in another email.

Sherry


> Anyway, is there a way to output the mat before it fails? We have been
> struggling to do that. We have TSSolve->SNES->FDColoringJacobian->A few
> levels of fieldsplit->failed Subblock matrix, which we want to analyze.
> (Sometimes it even happens in the second Newton iteration as the first one
> works okay.)
>
> Qi
>
>
>
> On Jun 10, 2022, at 8:11 AM, Matthew Knepley <knepley at gmail.com> wrote:
>
> 
> On Thu, Jun 9, 2022 at 5:20 PM Jorti, Zakariae via petsc-users <
> petsc-users at mcs.anl.gov> wrote:
>
>> Hi,
>>
>> I am solving non-linear problem that has 5 unknowns {ni, T, E, B, V}, and
>> for the preconditioning part, I am using a FieldSplit preconditioner. At
>> the last fieldsplit/level, we are left with a {B,V} block that tried to
>> precondition in 2 different ways:
>> a) SuperLU:
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_ksp_type preonly
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_pc_type lu
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_pc_factor_mat_solver_type
>> superlu_dist
>> b) a Schur-based fieldsplit preconditioner that uses SuperLU for both V
>> and B blocks:
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_ksp_type gmres
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_pc_type fieldsplit
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_pc_fieldsplit_type schur
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_pc_fieldsplit_schur_precondition
>> selfp -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_B_ksp_type
>> preonly -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_B_pc_type
>> lu
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_B_pc_factor_mat_solver_type
>> superlu_dist
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_V_ksp_type preonly
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_V_pc_type lu
>> -fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_fieldsplit_V_pc_factor_mat_solver_type
>> superlu_dist
>>
>> Option a) yields the following error:
>> "     Linear fieldsplit_ni_ solve converged due to CONVERGED_ATOL
>> iterations 0
>>         Linear fieldsplit_TEBV_fieldsplit_tau_ solve converged due to
>> CONVERGED_RTOL iterations 1
>>           Linear fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_EP_ solve
>> converged due to CONVERGED_RTOL iterations 5
>>           Linear fieldsplit_TEBV_fieldsplit_EBV_fieldsplit_BV_ solve did
>> not converge due to DIVERGED_PC_FAILED iterations 0
>>                          PC failed due to FACTOR_NUMERIC_ZEROPIVOT "
>> whereas options b) seems to be working well.
>> Is it possible that the SuperLU on the {V,B} block uses a reordering that
>> introduces a zero pivot or could there be another explanation for this
>> error?
>>
>
> I can at least come up with a case where this is true. Suppose you have
>
>   / A 0 \
>   \ 0  I /
>
> where A is rank deficient, but has a positive diagonal. SuperLU will fail
> since it is actually singular. However, your Schur complement might work
> since you use
> 'selfp' for the Schur preconditioner, and it just extracts the diagonal.
>
>   Thanks,
>
>      Matt
>
>
>> Many thanks.
>> Best,
>>
>> Zakariae
>>
>>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://www.cse.buffalo.edu/~knepley/
> <https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!HXCxUKc!3sg7TZfA3Z8m6foZ2wvm-LDZ6jI9C-kp1FVdsCbdXrDj--rBHyrWkz0akiZApgTAMNYce0yg5uta7w$>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220610/75a18aeb/attachment-0001.html>


More information about the petsc-users mailing list