[petsc-users] Improving efficiency of slepc usage

Jose E. Roman jroman at dsic.upv.es
Fri Aug 20 11:20:24 CDT 2021


Maybe too much fill-in during factorization. Try using an external linear solver such as MUMPS as explained in section 3.4.1 of SLEPc's users manual.

Jose


> El 20 ago 2021, a las 16:12, Matthew Knepley <knepley at gmail.com> escribió:
> 
> On Fri, Aug 20, 2021 at 6:55 AM dazza simplythebest <sayosale at hotmail.com> wrote:
> Dear Jose,
>     Many thanks for your response, I have been investigating this issue with a few more calculations 
> today, hence the slightly delayed response.
> 
> The problem is actually derived from a fluid dynamics problem, so to allow an easier exploration of things 
> I first downsized the resolution of the underlying fluid solver while keeping all the physical parameters
>  the same - i.e. I would get a smaller matrix that should be solving the same physical problem as the original
>  larger matrix but to lower accuracy.  
> 
> Results
> 
> Small matrix (N= 21168) - everything good!
> This converged when using the -eps_largest_real approach (taking 92 iterations for nev=10,
> tol= 5.0000E-06 and ncv = 300), and also when using the shift-invert approach, converging 
> very impressively in a single iteration ! Interestingly it did this both for a non-zero  -eps_target
>  and also for a zero  -eps_target.
> 
> Large matrix (N=50400)- works for -eps_largest_real , fails for st_type sinvert 
> I have just double checked again that the code does run properly when we use the -eps_largest_real 
> option - indeed I ran it with a small nev and large tolerance (nev = 4, tol= -eps_tol 5.0e-4 , ncv = 300)
> and with these parameters convergence was obtained in 164 iterations, which took 6 hours on the 
> machine I was running it on. Furthermore the eigenvalues seem to be ballpark correct; for this large
> higher resolution case (although with lower slepc tolerance) we obtain 1789.56816314173 -4724.51319554773i
>  as the eigenvalue with largest real part, while the smaller matrix (same physical problem but at lower resolution case)
> found this eigenvalue to be 1831.11845726501 -4787.54519511345i , which means the agreement is in line
> with expectations.
> 
> Unfortunately though the code does still crash though when I try to do shift-invert for the large matrix case ,
>  whether or not I use a non-zero  -eps_target. For reference this is the command line used :
> -eps_nev 10    -eps_ncv 300  -log_view -eps_view   -eps_target 0.1 -st_type sinvert -eps_monitor :monitor_output05.txt  
> To be precise the code crashes soon after calling EPSSolve (it successfully calls 
>  MatCreateVecs, EPSCreate,  EPSSetOperators, EPSSetProblemType and EPSSetFromOptions).
> By crashes I mean that I do not even get any error messages from slepc/PETSC, and do not even get the 
> 'EPS Object: 16 MPI processes' message - I simply get a  MPI/Fortran 'KILLED BY SIGNAL: 9 (Killed)' message
>  as soon as EPSsolve is called.
> 
> Hi Dan,
> 
> It would help track this error down if we had a stack trace. You can get a stack trace from the debugger. You run with
> 
>   -start_in_debugger
> 
> which should launch the debugger (usually), and then type
> 
>   cont
> 
> to continue, and then
> 
>   where
> 
> to get the stack trace when it crashes, or 'bt' on lldb.
> 
>   Thanks,
> 
>      Matt
>  
> Do you have any ideas as to why this larger matrix case should fail when using shift-invert but succeed when using 
> -eps_largest_real ? The fact that the program works and produces correct results 
> when using the -eps_largest_real  option suggests that there is probably nothing wrong with the specification 
> of the problem or the matrices ? It is strange how there is no error message from slepc / Petsc ... the 
> only idea I have at the moment is that perhaps max memory has been exceeded, which could cause such a sudden 
> shutdown? For your reference when running the large matrix case with the -eps_largest_real option I am using 
> about 36 GB of the 148GB available on this machine  - does the shift invert approach require substantially 
> more memory for example ?
> 
>   I would be very grateful if you have any suggestions to resolve this issue or even ways to clarify it further,
>  the performance I have seen with the shift-invert for the small matrix is so impressive it would be great to
>  get that working for the full-size problem.
> 
>    Many thanks and best wishes,
>                                   Dan.
> 
> 
> 
> From: Jose E. Roman <jroman at dsic.upv.es>
> Sent: Thursday, August 19, 2021 7:58 AM
> To: dazza simplythebest <sayosale at hotmail.com>
> Cc: PETSc <petsc-users at mcs.anl.gov>
> Subject: Re: [petsc-users] Improving efficiency of slepc usage
>  
> In A) convergence may be slow, especially if the wanted eigenvalues have small magnitude. I would not say 600 iterations is a lot, you probably need many more. In most cases, approach B) is better because it improves convergence of eigenvalues close to the target, but it requires prior knowledge of your spectrum distribution in order to choose an appropriate target.
> 
> In B) what do you mean that it crashes. If you get an error about factorization, it means that your A-matrix is singular, In that case, try using a nonzero target -eps_target 0.1
> 
> Jose
> 
> 
> > El 19 ago 2021, a las 7:12, dazza simplythebest <sayosale at hotmail.com> escribió:
> > 
> > Dear All,
> >             I am planning on using slepc to do a large number of eigenvalue calculations
> >  of a generalized eigenvalue problem, called from a program written in fortran using MPI.
> >  Thus far I have successfully installed the slepc/PETSc software, both locally and on a cluster,
> >  and on smaller test problems everything is working well; the matrices are efficiently and 
> > correctly constructed and slepc returns the correct spectrum. I am just now starting to move
> > towards now solving the full-size 'production run' problems, and would appreciate some 
> > general advice on how to improve the solver's performance.
> > 
> > In particular, I am currently trying to solve the problem Ax = lambda Bx whose matrices 
> > are of size 50000 (this is the smallest 'production run' problem I will be tackling), and are 
> > complex, non-Hermitian.  In most cases I aim to find the eigenvalues with the largest real part, 
> > although in other cases I will also be interested in finding the eigenvalues whose real part 
> > is close to zero.
> > 
> > A)
> > Calling slepc 's EPS solver with the following options:
> > 
> > -eps_nev 10   -log_view -eps_view -eps_max_it 600 -eps_ncv 140  -eps_tol 5.0e-6  -eps_largest_real -eps_monitor :monitor_output.txt
> > 
> > 
> > led to the code successfully running, but failing to find any eigenvalues within the maximum 600 iterations 
> > (examining the monitor output it did appear to be very slowly approaching convergence).
> > 
> > B)
> > On the same problem I have also tried a shift-invert transformation using the options
> > 
> > -eps_nev 10    -eps_ncv 140    -eps_target 0.0+0.0i  -st_type sinvert
> > 
> > -in this case the code crashed at the point it tried to call slepc, so perhaps I have incorrectly specified these options ?
> > 
> > 
> > Does anyone have any suggestions as to how to improve this performance ( or find out more about the problem) ?
> > In the case of A) I can see from watching the slepc   videos that increasing ncv 
> > may help, but I am wondering , since 600 is a large number of iterations, whether there 
> > maybe something else going on - e.g. perhaps some alternative preconditioner may help ?
> > In the case of B), I guess there must be some mistake in these command line options?
> >  Again, any advice will be greatly appreciated.
> >      Best wishes,  Dan.
> 
> 
> 
> -- 
> 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/



More information about the petsc-users mailing list