[Trilinos-Users] Large Generalized Symmetric Complex Eigenvalue Problem: From ARPACK to Anasazi?

Heroux, Michael A maherou at sandia.gov
Sun Aug 8 20:46:15 MDT 2010


Klaus,

Sorry about the problems you are having.  These errors should be fairly straight-forward to fix.  We will get to them as soon as possible and let you know when they are done.

Best regards,

Mike


On 8/7/10 11:11 AM, "Klaus Zimmermann" <klaus.zimmermann at physik.uni-freiburg.de> wrote:

Hello Heidi,

I followed our advise and started to implement the suggested software stack.
Unfortunately Tifpack seems to have some problems with complex numbers.
At various points in the code it tries to use operator< and/or operator> on them.
For example when I try to use ILUT I get [1]. The same problem exists with the other solvers,
which becomes apparent via use of the Tifpack::Factory [2].
I'd be happy with ILUT only for starters.

What can I do?

Thanks in advance,
Klaus

[1]:
------------------------------------
[kz38 at bwui helium1d_trilinos]$ mpic++ -I$HOME/usr/include main.cpp mat_he1d_floquet.cpp ~/usr/lib/libteuchos.a ~/usr/lib/libtpetra.a ~/usr/lib/libkokkos.a ~/usr/lib/libkokkosnodeapi.a ~/usr/lib/libteuchos.a ~/usr/lib/libtpi.a
/home/kz38/usr/include/Tifpack_ILUT_def.hpp: In member function ‘void Tifpack::ILUT<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:124:   instantiated from here
/home/kz38/usr/include/Tifpack_ILUT_def.hpp:363: error: no match for ‘operator<’ in ‘v < 0.0’
/home/kz38/usr/include/Tifpack_Heap.hpp: In member function ‘bool Tifpack::greater_indirect<Scalar, Ordinal>::operator()(const Ordinal&, const Ordinal&) const [with Scalar = std::complex<double>, Ordinal = int]’:
/opt/bwgrid/compiler/gnu/4.3.3/lib/gcc/x86_64-unknown-linux-gnu/4.3.3/../../../../include/c++/4.3.3/bits/stl_heap.h:180:  instantiated from ‘void std::__push_heap(_RandomAccessIterator, _Distance, _Distance, _Tp, _Compare) [with _RandomAccessIterator = __gnu_cxx::__normal_iterator<int*, std::vector<int, std::allocator<int> > >, _Distance = long int, _Tp = int, _Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/opt/bwgrid/compiler/gnu/4.3.3/lib/gcc/x86_64-unknown-linux-gnu/4.3.3/../../../../include/c++/4.3.3/bits/stl_heap.h:218:  instantiated from ‘void std::push_heap(_RAIter, _RAIter, _Compare) [with _RAIter = __gnu_cxx::__normal_iterator<int*, std::vector<int, std::allocator<int> > >, _Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/home/kz38/usr/include/Tifpack_Heap.hpp:72:   instantiated from ‘void Tifpack::add_to_heap(const Ordinal&, Teuchos::Array<T>&, SizeType&, Compare) [with Ordinal = int, SizeType = Tifpack::ILUT::compute::Tsize_t, Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/home/kz38/usr/include/Tifpack_ILUT_def.hpp:384:   instantiated from ‘void Tifpack::ILUT<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’
main.cpp:124:   instantiated from here
/home/kz38/usr/include/Tifpack_Heap.hpp:45: error: no match for ‘operator>’ in ‘((const Teuchos::Array<std::complex<double> >*)((const Tifpack::greater_indirect<std::complex<double>, int>*)this)->Tifpack::greater_indirect<std::complex<double>, int>::m_vals)->Teuchos::Array<T>::operator[] [with T = std::complex<double>](((long int)lhs)) > ((const Teuchos::Array<std::complex<double> >*)((const Tifpack::greater_indirect<std::complex<double>, int>*)this)->Tifpack::greater_indirect<std::complex<double>, int>::m_vals)->Teuchos::Array<T>::operator[] [with T = std::complex<double>](((long int)rhs))’
------------------------------------

[2]:
------------------------------------
[kz38 at bwui helium1d_trilinos]$ mpic++ -I$HOME/usr/include main.cpp mat_he1d_floquet.cpp ~/usr/lib/libteuchos.a ~/usr/lib/libtpetra.a ~/usr/lib/libkokkos.a ~/usr/lib/libkokkosnodeapi.a ~/usr/lib/libteuchos.a ~/usr/lib/libtpi.a
/home/kz38/usr/include/Tifpack_Chebyshev_def.hpp: In member function ‘void Tifpack::Chebyshev<MatrixType>::describe(Teuchos::FancyOStream&, Teuchos::EVerbosityLevel) const [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_Chebyshev_def.hpp:540: error: no match for ‘operator>’ in ‘myMinVal > DiagView.Teuchos::ArrayRCP<T>::operator[] [with T = const std::complex<double>](i)’
/home/kz38/usr/include/Tifpack_Chebyshev_def.hpp:541: error: no match for ‘operator<’ in ‘myMaxVal < DiagView.Teuchos::ArrayRCP<T>::operator[] [with T = const std::complex<double>](i)’
/home/kz38/usr/include/Tifpack_Chebyshev_def.hpp: In member function ‘void Tifpack::Chebyshev<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_Chebyshev_def.hpp:408: error: no match for ‘operator<’ in ‘Teuchos::ScalarTraits<std::complex<_Tp> >::magnitude [with T = double](diag) < ((Tifpack::Chebyshev<Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> > >*)this)->Tifpack::Chebyshev<Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> > >::MinDiagonalValue_’
/home/kz38/usr/include/Tifpack_Relaxation_def.hpp: In member function ‘void Tifpack::Relaxation<MatrixType>::describe(Teuchos::FancyOStream&, Teuchos::EVerbosityLevel) const [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_Relaxation_def.hpp:828: error: no match for ‘operator>’ in ‘myMinVal > DiagView.Teuchos::ArrayRCP<T>::operator[] [with T = std::complex<double>](i)’
/home/kz38/usr/include/Tifpack_Relaxation_def.hpp:829: error: no match for ‘operator<’ in ‘myMaxVal < DiagView.Teuchos::ArrayRCP<T>::operator[] [with T = std::complex<double>](i)’
/home/kz38/usr/include/Tifpack_Relaxation_def.hpp: In member function ‘void Tifpack::Relaxation<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_Relaxation_def.hpp:363: error: no match for ‘operator<’ in ‘Teuchos::ScalarTraits<std::complex<_Tp> >::magnitude [with T = double](diag) < ((Tifpack::Relaxation<Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> > >*)this)->Tifpack::Relaxation<Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> > >::MinDiagonalValue_’
/home/kz38/usr/include/Tifpack_RILUK_def.hpp: In member function ‘void Tifpack::RILUK<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_RILUK_def.hpp:364: error: no match for ‘operator>’ in ‘Teuchos::ScalarTraits<std::complex<_Tp> >::magnitude [with T = double](DV.Teuchos::ArrayRCP<T>::operator[] [with T = std::complex<double>](((long int)i))) > MaxDiagonalValue’
/home/kz38/usr/include/Tifpack_RILUK_def.hpp:365: error: no match for ‘operator<’ in ‘DV.Teuchos::ArrayRCP<T>::operator[] [with T = std::complex<double>](((long int)i)) < 0’
/home/kz38/usr/include/Tifpack_ILUT_def.hpp: In member function ‘void Tifpack::ILUT<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_ILUT_def.hpp:363: error: no match for ‘operator<’ in ‘v < 0.0’
/home/kz38/usr/include/Tifpack_RILUK_def.hpp: In member function ‘void Tifpack::RILUK<MatrixType>::initAllValues(const Tpetra::RowMatrix<typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type>&) [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’:
/home/kz38/usr/include/Tifpack_RILUK_def.hpp:144:   instantiated from ‘void Tifpack::RILUK<MatrixType>::initialize() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_RILUK_def.hpp:200: error: no match for ‘operator<’ in ‘InV.Teuchos::Array<T>::operator[] [with T = std::complex<double>](((long int)j)) < 0.0’
/home/kz38/usr/include/Tifpack_Heap.hpp: In member function ‘bool Tifpack::greater_indirect<Scalar, Ordinal>::operator()(const Ordinal&, const Ordinal&) const [with Scalar = std::complex<double>, Ordinal = int]’:
/opt/bwgrid/compiler/gnu/4.3.3/lib/gcc/x86_64-unknown-linux-gnu/4.3.3/../../../../include/c++/4.3.3/bits/stl_heap.h:180:  instantiated from ‘void std::__push_heap(_RandomAccessIterator, _Distance, _Distance, _Tp, _Compare) [with _RandomAccessIterator = __gnu_cxx::__normal_iterator<int*, std::vector<int, std::allocator<int> > >, _Distance = long int, _Tp = int, _Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/opt/bwgrid/compiler/gnu/4.3.3/lib/gcc/x86_64-unknown-linux-gnu/4.3.3/../../../../include/c++/4.3.3/bits/stl_heap.h:218:  instantiated from ‘void std::push_heap(_RAIter, _RAIter, _Compare) [with _RAIter = __gnu_cxx::__normal_iterator<int*, std::vector<int, std::allocator<int> > >, _Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/home/kz38/usr/include/Tifpack_Heap.hpp:72:   instantiated from ‘void Tifpack::add_to_heap(const Ordinal&, Teuchos::Array<T>&, SizeType&, Compare) [with Ordinal = int, SizeType = Tifpack::ILUT::compute::Tsize_t, Compare = Tifpack::greater_indirect<std::complex<double>, int>]’
/home/kz38/usr/include/Tifpack_ILUT_def.hpp:384:   instantiated from ‘void Tifpack::ILUT<MatrixType>::compute() [with MatrixType = Tpetra::CrsMatrix<std::complex<double>, int, int, Kokkos::TPINode, Kokkos::DefaultSparseOps<std::complex<double>, int, Kokkos::TPINode> >]’
main.cpp:129:   instantiated from here
/home/kz38/usr/include/Tifpack_Heap.hpp:45: error: no match for ‘operator>’ in ‘((const Teuchos::Array<std::complex<double> >*)((const Tifpack::greater_indirect<std::complex<double>, int>*)this)->Tifpack::greater_indirect<std::complex<double>, int>::m_vals)->Teuchos::Array<T>::operator[] [with T = std::complex<double>](((long int)lhs)) > ((const Teuchos::Array<std::complex<double> >*)((const Tifpack::greater_indirect<std::complex<double>, int>*)this)->Tifpack::greater_indirect<std::complex<double>, int>::m_vals)->Teuchos::Array<T>::operator[] [with T = std::complex<double>](((long int)rhs))’
------------------------------------

On 07.08.10 00:23, Thornquist, Heidi K wrote:
Re: [Trilinos-Users] Large Generalized Symmetric Complex Eigenvalue Problem: From ARPACK to Anasazi? Hi Klaus,

For now we’d suggest Anasazi (BKS), Tpetra, Belos (GMRES), Tifpack (for preconditioning).  Since Tifpack is a
domain decomposition incomplete LU, it won’t take nearly as much memory as Pardiso.  This is because it
is not performing an full factorization of the matrix, just an incomplete LU factorization of the parts of the domain that
are on each processor.  When you are using ARPACK, I’m assuming that you are using a spectral transformation,
which is why you need the linear solver.  If so, which linear solver did you use with ARPACK?

Thanks,
Heidi


On 8/6/10 4:13 PM, "Klaus Zimmermann" <klaus.zimmermann at physik.uni-freiburg.de> <mailto:klaus.zimmermann at physik.uni-freiburg.de>  wrote:


Mike,

thanks for your quick reply!
I am glad to hear of these improvements. If I understand you correctly
for now you'd suggest
Anasazi, Tpetra, Belos, right?
Then there is the question of the preconditioner: In our legacy code we
did a LU factorization
(via pardiso) which quickly led to storage requirements on the order of
50 gb. This is unfortunately
prohibitive for the 16 gb of ram/node on the grid. We used to use out of
core functionality for the
solver.
Is such out of core stuff available for Belos as well, or will Tifpack
reduce the amount of memory
needed so much I don't have to worry?
For reference  the matrices we deal with right now are on the order of
n=1e8, complex symmetric
with about 0.5% non zero elements.

Klaus

On 06.08.10 23:53, Heroux, Michael A wrote:
> Klaus,
>
> The software stack that supports complex is growing quickly.  Belos already
> supports complex for Krylov methods for linear systems.  There is a new
> package called Tifpack that supports complex to do domain decomposition with
> ILU and related preconditioners.  So, using Trilinos 10.4 (the current
> release), there is a better answer than Komplex for iterative linear
> solvers.
>
> In the next version of Trilinos (Oct 2010), we plan to have Amesos2 and KLU
> 2, which support complex.  How much of the linking support to complex
> versions of MUMPS, UMFPACK, etc. will be done by October is not clear right
> now, but we are moving in the right direction.
>
> BTW, all of the support for complex comes from using templated scalar types
>
> Mike
>
> On 8/6/10 4:41 PM, "Klaus Zimmermann"
> <klaus.zimmermann at physik.uni-freiburg.de> <mailto:klaus.zimmermann at physik.uni-freiburg.de>   wrote:
>
>
>> Hi All,
>>
>> for my ab initio quantum calculations I have to find some eigenvalues
>> for large, sparse complex symmetric generalized eigenvalue problems.
>> I would really like to move from our legacy ARPACK fortran code to a
>> more modern code base
>> in order to make use of our local MPI based grid. I am a bit overwhelmed
>> by the choices. But I now think
>> I boiled my possibilities down to Anasazi or SLEPc. I'd prefer to go the
>> C++ way.
>> However I understand the support for complex systems in the linear
>> solver part is limited.
>>
>> What would be the best way to go about tackling this kind of problem?
>> Right now I am thinking: Anasazi<Tpetra, BlockKrylovSchur>  +
>> (Komplex+Amesos+MUMPS).
>> Is this reasonable?
>>
>> Thanks for any hints in advance!
>>
>> Regards,
>> Klaus Zimmermann
>>
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at software.sandia.gov
>> http://software.sandia.gov/mailman/listinfo/trilinos-users
>>
>
>


_______________________________________________
Trilinos-Users mailing list
Trilinos-Users at software.sandia.gov
 http://software.sandia.gov/mailman/listinfo/trilinos-users






-------------- next part --------------
An HTML attachment was scrubbed...
URL: https://software.sandia.gov/pipermail/trilinos-users/attachments/20100808/b27494c1/attachment-0001.html 


More information about the Trilinos-Users mailing list