[petsc-users] KrylovSchur solver diverges when setting EPS_GHEP
    Denis Davydov 
    davydden at gmail.com
       
    Tue Nov 25 16:00:25 CST 2014
    
    
  
> 
> sinvert should give you more that 2 eigenvalues if you set eps_nev > 2. You may need to increase the number of iterations with eps_max_it
i would expect that it should, but it is not the case. I tried setting number of iterations 10 times the number of DoFs, and no change.
It does not seem to be related. Here is an output without shift-and-invert
EPS Object: 1 MPI processes
  type: krylovschur
    Krylov-Schur: 50% of basis vectors kept after restart
  problem type: generalized symmetric eigenvalue problem
  selected portion of the spectrum: smallest real parts
  computing true residuals explicitly
  number of eigenvalues (nev): 5
  number of column vectors (ncv): 20
  maximum dimension of projected problem (mpd): 20
  maximum number of iterations: 2890
  tolerance: 1e-07
  convergence test: absolute
BV Object: 1 MPI processes
  type: svec
  21 columns of global length 289
  orthogonalization method: classical Gram-Schmidt
  orthogonalization refinement: if needed (eta: 0.7071)
  non-standard inner product
  Mat Object:   1 MPI processes
    type: seqaij
    rows=289, cols=289
    total: nonzeros=1913, allocated nonzeros=5491
    total number of mallocs used during MatSetValues calls =0
      not using I-node routines
DS Object: 1 MPI processes
  type: hep
  solving the problem with: Implicit QR method (_steqr)
ST Object: 1 MPI processes
  type: shift
  shift: 0
  number of matrices: 2
  all matrices have different nonzero pattern
  KSP Object:  (st_)   1 MPI processes
    type: preonly
    maximum iterations=10000, initial guess is zero
    tolerances:  relative=1e-08, absolute=1e-50, divergence=10000
    left preconditioning
    using NONE norm type for convergence test
  PC Object:  (st_)   1 MPI processes
    type: redundant
      Redundant preconditioner: First (color=0) of 1 PCs follows
      KSP Object:      (st_redundant_)       1 MPI processes
        type: preonly
        maximum iterations=10000, initial guess is zero
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
        left preconditioning
        using NONE norm type for convergence test
      PC Object:      (st_redundant_)       1 MPI processes
        type: lu
          LU: out-of-place factorization
          tolerance for zero pivot 2.22045e-14
          matrix ordering: nd
          factor fill ratio given 5, needed 2.99686
            Factored matrix follows:
              Mat Object:               1 MPI processes
                type: seqaij
                rows=289, cols=289
                package used to perform factorization: petsc
                total: nonzeros=5733, allocated nonzeros=5733
                total number of mallocs used during MatSetValues calls =0
                  not using I-node routines
        linear system matrix = precond matrix:
        Mat Object:         1 MPI processes
          type: seqaij
          rows=289, cols=289
          total: nonzeros=1913, allocated nonzeros=5491
          total number of mallocs used during MatSetValues calls =0
            not using I-node routines
    linear system matrix = precond matrix:
    Mat Object:     1 MPI processes
      type: seqaij
      rows=289, cols=289
      total: nonzeros=1913, allocated nonzeros=5491
      total number of mallocs used during MatSetValues calls =0
        not using I-node routines
> 
>> 
>> 2) Another peculiarity is that `-eps_smallest_real` and `-eps_target 0 -st_type sinvert` return different sets of eigenvalues. 
> 
> The matrix from the 2D Laplace operator is positive definite, so both should give the same eigenvalues.
> 
>> In the latter case there are degenerate eigenvalues. 
>> Those are consistent with the results given by ARPACK with shift and invert around 0.
>> The matrices I sent you originally correspond to the eigenvalues of the 2D Laplace operator on a uniform mesh with 256 cells.
>> If more refined mesh is used (1024 cells instead of 256), same set with degenerate eigenvalues is returned in both cases.
>> This is not directly related to SLEPc and is a question out of curiosity.
> 
> Some eigenvalues of the 2D Laplacian have multiplicity 2. The default SLEPc solver may not return the two copies of the eigenvalue, because the second copy appears much later and the iteration stops as soon as nev eigenvalues have been computed.
i see… 
so technically, one could try asking for more eigenvectors to resolve those degenerate eigenvalues?
> 
>> 
>> 3) It is not stated in the documentation explicitly, but I suppose the residual discussed in ‘SolverControl’ section 
>> always corresponds to the direct problem (even in case when, say, shift-and-invert is applied) and so 
>> EPSComputeRelativeError and EPSComputeResidualNorm?
> 
> The convergence criterion is applied to the transformed problem (shift-and-invert).
understood, thanks.
Regards,
Denis.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20141125/346a24c6/attachment-0001.html>
    
    
More information about the petsc-users
mailing list