[Trilinos-Users] Applying Anasazi::BlockKrylovSchurSolMgr

Rutherford, Joseph M jmruther at illinois.edu
Mon Dec 14 09:49:28 EST 2015


Alicia et al,

Thanks for your very helpful reply.  I extended the size of my 4x4 system to 10x10 by increasing the dimension size in the serial Tpetra::Map instance I’m using for both domain and range. So now I still have a pair of rank 4 matrices, and the A x = \lambda B x system still should have the same eigenvalues and the eigenvectors should just have zero padding.

The good news is that your suggestion worked very well; the system is now satisfied to move into the solve stage. The bad news is that the solver seems to be getting stuck in a loop in the first iteration.  Since I can SVD a dense random 100x100 system nearly instantaneously with NumPy, I think I must still be using BKS incorrectly if its first iteration for 10x10 is still underway after several minutes of full CPU usage.  Do you have any suggestions for how to dig further into this?  Debugging indicates I’m stuck in the while(doGramSchmidt){…} loop of int SVQBOrthoManager<ScalarType, MV, OP>::findBasis(…).  Terminal output below.

Thanks again for getting me past that last hurdle.

Joe

<terminal>
$ ./bin/test_eigensolver.exe
Running 1 test case...
Debugging checks: iteration 0: after initialize()
>> Error in F^H M F == I  : 1.11e+01

================================================================================

                         BlockKrylovSchur Solver Status

The solver is initialized.
The number of iterations performed is 0
The block size is         1
The number of blocks is   5
The current basis size is 0
The number of auxiliary vectors is 0
The number of operations Op*x   is 0

CURRENT RITZ VALUES
   [ NONE COMPUTED ]

================================================================================

- StatusTestOutput: Failed
  (Num calls,Mod test,State test): (1, 1, Passed Failed Undefined )
   - StatusTestCombo: Failed
     - StatusTestWithOrdering: Undefined
       Quorum: 4
       Auxiliary values: [empty]
       - StatusTestResNorm: Failed
         (Tolerance,WhichNorm,Scaled,Quorum): (1.000000e-03,RITZRES_2NORM,true,4)
         Which vectors: [empty]
</terminal>


From: Alicia Klinvex [mailto:aklinvex at purdue.edu]
Sent: Sunday, December 13, 2015 12:33 AM
To: Rutherford, Joseph M <jmruther at illinois.edu>
Cc: trilinos-users at trilinos.org
Subject: Re: [Trilinos-Users] Applying Anasazi::BlockKrylovSchurSolMgr

Hello Joseph,

BKS tries to build a subspace of dimension block_size * num_blocks.  If the dimension of your problem n is less than that subspace dimension, it's impossible to orthogonalize the basis vectors (since you can't generate 10 orthogonal vectors within a subspace of dimension 4).  Perhaps we should modify the error message to make this clearer.

Can you try a slightly larger example, perhaps 10x10, with block_size 1 and num_blocks < 10?

Best wishes,
Alicia

On Sat, Dec 12, 2015 at 1:57 PM, Rutherford, Joseph M <jmruther at illinois.edu<mailto:jmruther at illinois.edu>> wrote:
I’m having some difficulty using Anasazi to compute the smallest eigenvalues of a (non-Hermitian) complex-valued linear system.  Before addressing a general case,  I have written a simple test using two Tpetra::CrsMatrix<std::complex<float>,int,int> instances with reference data – this is a very small 4x4 dense case with reference data available.  (I can share the example code if needed.)

Unfortunately, I am having a hard time getting the SolMgr to drive the solution.  As I’m solving a reference 4x4 system with 4 nondegenerate eigenvalues, I’m fairly confident I don’t have a pathological case, although it is much smaller than the usual Anasazi use case.  Could you please look at the error message reported to my console (below)?  I’ve not found a clearly-aligned example in the Anasazi distribution, so I’m hopeful this is just a trivial case that is easily fixed.

Sincerely yours,

Joe Rutherford

<console>
$ ./bin/test_eigensolver.exe
Running 1 test case...
unknown location(0): fatal error: in "eigensolver_dense_4x4_complex": std::invalid_argument: /home/jmruther/branches/bean/externals/trilinos/include/AnasaziBlockKrylovSchurSolMgr.hpp:344:

Throw number = 2

Throw test that evaluated to true: static_cast<ptrdiff_t>(_numBlocks)*_blockSize > MVT::GetGlobalLength(*_problem->getInitVec())

Anasazi::BlockKrylovSchurSolMgr: Potentially impossible orthogonality requests. Reduce basis size.
/home/jmruther/branches/bean/source/tests/eigensolver.cpp(57): last checkpoint: "eigensolver_dense_4x4_complex" entry.

*** 1 failure is detected in the test module "Eigensolver"
$
</console>



_______________________________________________
Trilinos-Users mailing list
Trilinos-Users at trilinos.org<mailto:Trilinos-Users at trilinos.org>
https://trilinos.org/mailman/listinfo/trilinos-users<https://urldefense.proofpoint.com/v2/url?u=https-3A__trilinos.org_mailman_listinfo_trilinos-2Dusers&d=BQMFaQ&c=8hUWFZcy2Z-Za5rBPlktOQ&r=Ng9lEtA9hn_MbIhqAtV5991RoWo8JGET5un-J3-Ub7o&m=a92OB0Je5c97qWnvqkK3eKwfHg65YXj6gsz4PO4b8P0&s=-Aqn6cS4l23CGmblqMY86E-cVM9F7hK4MpfCFtvnIkA&e=>

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20151214/c55ab0f7/attachment-0001.html>


More information about the Trilinos-Users mailing list