Run your code with -eps_view_mat0 binary:amatrix.bin
This will save the A matrix of the EPS solver in file "amatrix.bin" at the end 
of EPSSolve().
I cannot say which is the best method, it depends on the distribution of the 
spectrum. Let's see how it goes.

Jose


> El 11 oct 2018, a las 20:13, Moritz Cygorek <mcygo...@uottawa.ca> escribió:
> 
> Thank you for the response. This example is exactly the kind of thing I 
> thought about.
> But, as you say, convergence is indeed not better. 
> 
> 
> What's the best format to send the matrix to slepc-maint? 
> My small test case is a sparse matrix of dimension 30720. I currently read it 
> from three files (containing IA, JA, and the values in CSR format). Because 
> everything is contained in a larger framework, it's not very easy to extract 
> the reader.  Is there a generic way to read in sparse matrices from files for 
> PETSc?
> 
> 
> 
> Maybe you can already give me a hint about the best method if I tell you the 
> properties of the spectrum:
> The problem is the calculation of electronic states in a quantum dot.  There 
> is a band gap from around 0 to 2 and I am interested in the first 20 or 40 
> eigenvectors above and below the gap.
> The first few states are separated by about 0.1 but for higher states the 
> energies come closer and closer together.
> Additionally, there are a number of states with energies around 1000, which 
> are artificial and originate from the way we treat the boundary conditions. 
> Also, I know that all states come in pairs with the same energy (Kramers 
> degeneracy).
> 
> For the small test case, the conduction band states, i.e. eigenvectors with 
> energies close to 2, converge very fast (about 5 min on a laptop computer). 
> However, the states with energies around 0 converge much more slowly and 
> that's one of my major problems. 
> For those states harmonic extraction seems to be better suited, but I have 
> the impression that it is not extremely stable. For example, applied to the 
> states close to 2 I see that some states are skipped, which can be seen by 
> the fact that the degeneracies are sometimes wrong. 
> Also, with harmonic extraction, the program sometimes stops claiming the 
> number of requested eigenpairs are reached, but the calculated relative error 
> of most states is way larger than the tolerance.
> 
> Maybe you know from experience which method is better suited to tackle these 
> kinds of problems?
> 
> 
> Eventually, I intend to do calculates with dimensions of ~ 10 million 
> distributed of a few 100 CPUs.
> 
> Regards,
> Moritz
> 
> 
> 
> 
> From: Jose E. Roman <jro...@dsic.upv.es>
> Sent: Thursday, October 11, 2018 5:55 AM
> To: Matthew Knepley; Moritz Cygorek
> Cc: PETSc; Carmen Campos
> Subject: Re: [petsc-users] STFILTER in slepc
>  
> The filter technique must be used with an interval in EPS, e.g.
>   -eps_interval 2.,2.7
> (This interval will be passed to STFILTER, so no need to specify 
> -st_filter_interval).
> Therefore, it does not require -eps_target_magnitude or similar.
> 
> Regarding the simpler (A-tau*I)^2, we call it "spectrum folding" (see section 
> 3.4.6 of SLEPc's manual) and it is implemented in ex24.c
> http://slepc.upv.es/documentation/current/src/eps/examples/tutorials/ex24.c.html
> 
> My experience is that spectrum folding will not give good convergence except 
> for easy problems. For more difficult problems with clustered eigenvalues you 
> need a polynomial filter. The technique of polynomial filter trades 
> iterations of the Krylov subspace for matrix-vector products required by a 
> high-degree polynomial.
> 
> If you want, send us the matrix to slepc-maint and we will have a try.
> 
> Jose
>  
> 
> > El 11 oct 2018, a las 2:43, Matthew Knepley <knep...@gmail.com> escribió:
> > 
> > On Wed, Oct 10, 2018 at 8:41 PM Moritz Cygorek <mcygo...@uottawa.ca> wrote:
> > Thank you very much. Apparently, I've misunderstood what the filter 
> > actually does. I thought about the much simpler process, where you 
> > diagonalize
> > 
> > 
> > 
> > -(A- tau*I)^2 +offset*I
> > 
> > 
> > 
> > where tau is my target an offset is large enough so that the global maximum 
> > is reached for eigenvalues around tau.
> > 
> > 
> > Is this different from -eps_target_magnitude?
> > 
> >   Thanks,
> > 
> >     Matt
> >  
> > Then you look for the largest eigenvalue of the modified problem and either 
> > calculate the Ritz value of the original matrix or calculate back from the 
> > eigenvalues of the modified problem.
> > 
> > 
> > 
> > Now, it looks to me like -st_type filter activates something like the 
> > package FILTLAN.
> > 
> > 
> > 
> > I guess I can define a MatShell to do the thing I intended in the first 
> > place.
> > 
> > But, I guess, this is a common thing, so I am wondering whether it is 
> > already implemented somewhere and I just didn't find it in the 
> > documentation.  Can you say something about this?
> > 
> > 
> > 
> > Regards,
> > 
> > Moritz
> > 
> > 
> > 
> > 
> > 
> > From: Jose E. Roman <jro...@dsic.upv.es>
> > Sent: Wednesday, October 10, 2018 3:48 PM
> > To: Moritz Cygorek
> > Cc: petsc-users@mcs.anl.gov
> > Subject: Re: [petsc-users] STFILTER in slepc
> >  
> > This type of method requires a very high degree polynomial; suggest using 
> > degree=100 at least (this is the default value), but larger values may be 
> > necessary. Also, for this particular filter the "range" must be 
> > approximately equal to the numerical range; if you have no clue where your 
> > first and last eigenvalues are, you may use EPSSolve() calls with 
> > EPS_LARGEST_REAL and EPS_SMALLEST_REAL.
> > 
> > Jose
> > 
> > > El 10 oct 2018, a las 21:10, Moritz Cygorek <mcygo...@uottawa.ca> 
> > > escribió:
> > > 
> > > Thank you for the fast reply. 
> > > 
> > > I've tried running my program (using the defaul Krylov-Schur method for 
> > > sparse MPI matrices) with the additional options:
> > > 
> > > -st_type filter -st_filter_degree 2 -st_filter_interval 2.,2.7  
> > > -st_filter_range -2000,2000
> > > 
> > > and I get the following error message:
> > > 
> > > [0]PETSC ERROR: STFILTER cannot get the filter specified; please adjust 
> > > your filter parameters (e.g. increasing the polynomial degree)
> > > ....
> > > [0]PETSC ERROR: #1 FILTLAN_GetIntervals() line 451 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/sys/classes/st/impls/filter/filtlan.c
> > > [0]PETSC ERROR: #2 STFilter_FILTLAN_setFilter() line 1016 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/sys/classes/st/impls/filter/filtlan.c
> > > [0]PETSC ERROR: #3 STSetUp_Filter() line 42 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/sys/classes/st/impls/filter/filter.c
> > > [0]PETSC ERROR: #4 STSetUp() line 271 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/sys/classes/st/interface/stsolve.c
> > > [0]PETSC ERROR: #5 EPSSetUp() line 263 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/eps/interface/epssetup.c
> > > [0]PETSC ERROR: #6 EPSSolve() line 135 in 
> > > /home/applications/sources/libraries/slepc-3.9.2/src/eps/interface/epssolve.c
> > > 
> > > 
> > > 
> > > Do you have a clue what I've missed?
> > > 
> > > 
> > > Moritz
> > > 
> > > 
> > > From: Jose E. Roman <jro...@dsic.upv.es>
> > > Sent: Wednesday, October 10, 2018 2:30 PM
> > > To: Moritz Cygorek
> > > Cc: petsc-users@mcs.anl.gov
> > > Subject: Re: [petsc-users] STFILTER in slepc
> > >  
> > > 
> > > > El 10 oct 2018, a las 19:54, Moritz Cygorek <mcygo...@uottawa.ca> 
> > > > escribió:
> > > > 
> > > > Hi,
> > > > 
> > > > in the list of changes to SLEPc version 3.8, it is stated that there is 
> > > > a preliminary implementation of polynomial filtering using STFILTER. 
> > > > 
> > > > Because I am struggling to obtain interior eigenvalues and harmonic 
> > > > extraction seems not to be stable enough in my case, I wanted to give 
> > > > it a try, but I could not find any documentation yet.
> > > > 
> > > > Does anybody have an example of how to use STFILTER or any 
> > > > documentation about it?
> > > > 
> > > > Thanks in advance,
> > > > Moritz
> > > 
> > > There are no examples. You just set the type to STFILTER and set some 
> > > parameters such as the interval of interest or the polynomial degree. See 
> > > functions starting with STFilter 
> > > here:http://slepc.upv.es/documentation/current/docs/manualpages/ST/index.html
> > > 
> > > In some problems it works well, but don't expect too much. It is still in 
> > > our to-do list to make it more usable. It will be good to have your 
> > > feedback. If you want, send results to slepc-maint, maybe we can help 
> > > tuning the parameters.
> > > 
> > > Jose
> > 
> > 
> > 
> > -- 
> > 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/

Reply via email to