[petsc-users] How to compute the condition number of SchurComplementMat preconditioned with PCSHELL.

Barry Smith bsmith at petsc.dev
Thu Nov 11 09:14:26 CST 2021


  I think there is likely a problem somewhere along the line in the libraries or how you are using them that causes the problem. You can run with a small problem and explicit form the inner operators and use direct solvers for the inner inverses to see what happens to the convergence of the computation of the eigenvalues. You can also use -ksp_view_eigenvalues_explicit and things like -ksp_view_mat_explicit and -ksp_view_preconditioned_operator_explicit to dump representations of the operators (for small problems) and load them in Matlab to get the eigenvalues and see if they match what the PETSc eigenvalue computation gives.


> On Nov 11, 2021, at 4:39 AM, Vladislav Pimanov <Vladislav.Pimanov at skoltech.ru> wrote:
> 
> Dear Barry, 
> 
> Let S = B A^{-1} B^T be the Schur complement, and \hat{S} = B diag(A)^{-1} B^T denotes the preconditioner.
> I also tried KSPComputeExtremeSingularValues for rtol(A) = rtol(\hat{S}) = 1e-14 as you suggested. However, the results did not change much compared to rtol(A) = rtol(S) = rtol(\hat{S}). 
> 
> Additionally, I checked the identity preconditioner instead of \hat{S}. In this case, the method firstly seems to converge when reaching rtol(S)=1e-6, but then, with a further decrease of rtol(S), it begins to diverge slightly though.
> 
> I will try generalised EVP.
> 
> Thanks,
> Vlad
> 
> От: Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>>
> Отправлено: 10 ноября 2021 г. 17:45:55
> Кому: Vladislav Pimanov
> Копия: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>
> Тема: Re: [petsc-users] How to compute the condition number of SchurComplementMat preconditioned with PCSHELL.
>  
> 
>> P is a diffusion matrix, which itself is inverted by KSPCG)
> 
> This worries me. Unless solved to full precision the action of solving with CG is not a linear operator in the input variable b, this means that the action of your Schur complement is not a linear operator and so iterative eigenvalue algorithms may not work correctly (even if the linear system seems to be converging). 
> 
> I suggest you try a KSP tolerance of 10^-14 for the inner KSP solve when you attempt the eigenvalue computation.
> 
> Barry
> 
> 
>> On Nov 10, 2021, at 9:09 AM, Vladislav Pimanov <Vladislav.Pimanov at skoltech.ru <mailto:Vladislav.Pimanov at skoltech.ru>> wrote:
>> 
>> Great thanks, Matt!
>> The second option is what I was looking for.
>> 
>> Best Regards,
>> Vlad
>> От: Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>>
>> Отправлено: 10 ноября 2021 г. 16:45:00
>> Кому: Vladislav Pimanov
>> Копия: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>
>> Тема: Re: [petsc-users] How to compute the condition number of SchurComplementMat preconditioned with PCSHELL.
>>  
>> On Wed, Nov 10, 2021 at 8:42 AM Vladislav Pimanov <Vladislav.Pimanov at skoltech.ru <mailto:Vladislav.Pimanov at skoltech.ru>> wrote:
>> Dear PETSc community,
>> 
>> 
>> I wonder if you could give me a hint on how to compute the condition number of a preconditioned matrix in a proper way.
>> I have a MatSchurComplement matrix S and a preconditioner P of the type PCSHELL (P is a diffusion matrix, which itself is inverted by KSPCG).
>> I tried to compute the condition number of P^{-1}S "for free" during the outer PCG procedure using KSPComputeExtremeSingularValues() routine.
>> Unfortunately, \sigma_min does not converge even if the solution is computed with very high precision.
>> I also looked at SLEPc interface, but did not realised how PC should be included.
>> 
>> You can do this at least two ways:
>> 
>>   1) Make a MatShell for P^{-1} S. This is easy, but you will not be able to use any factorization-type PC on that matrix.
>> 
>>   2) Solve instead the generalized EVP, S x = \lambda P x. Since you already have P^{-1}, this should work well.
>> 
>>   Thanks,
>> 
>>      Matt
>> Thanks!
>> 
>> Sincerely,
>> Vladislav Pimanov
>> 
>> 
>> 
>> -- 
>> 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/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20211111/0f337228/attachment-0001.html>


More information about the petsc-users mailing list