[Trilinos-Users] Multiplying an Epetra_CrsMatrix by another Epetra_CrsMatrix

Heroux, Michael A maherou at sandia.gov
Tue May 25 09:27:52 MDT 2010


Mehmet,

Since B is diagonal, I would store B^{-1/2} as an Epetra_Vector, then use the RightScale and LeftScale methods available in the Epetra_CrsMatrix class:

http://trilinos.sandia.gov/packages/docs/r10.2/packages/epetra/doc/html/classEpetra__CrsMatrix.html#a73fa7db6dc1738bcb033d36362ffa36

http://trilinos.sandia.gov/packages/docs/r10.2/packages/epetra/doc/html/classEpetra__CrsMatrix.html#1d51803cf77bd9ee904952aaf03d8b8b

These methods will scale the matrix A in-place.

There is also a Multiply (or you probably want ReciprocalMultiply) method on the MultiVector class that will scale your solution back to the original space:

http://trilinos.sandia.gov/packages/docs/r10.2/packages/epetra/doc/html/classEpetra__MultiVector.html#7fc1add99930d33c11d34d198dab9656

Mike


On 5/25/10 3:34 AM, "Mehmet Salih YILDIRIM" <linux at isadamlari.org> wrote:

Hi all!

I have very large Epetra_CrsMatrix matrices and I have to multiply one with another. That is, I'm actually trying to transform a generalized eigenproblem into a standard one, and like most of you already know, my problem is
A * x = B * x * \lambda

should be transformed into something like:
B^{-1/2} * A * B{-1/2} * y = y * \lambda

It is easy to calculate B^{-1/2} since B is diagonal, but I couldn't figure out how to do B^{-1/2} * A * B^{-1/2}  multiplication because I am not sure if Epetra_CrsMatrix::Multiply satisfies my need here since its second and third arguements are Epetra_Vector or Epetra_MultiVector and they are not in the same class hierarchy with Epetra_CrsMatrix.

On the other hand, I have to keep B^{-1/2} in the memory for transofrming the solution back into generalized problem's solution, but i don't need to store A in the memory. A is a Teuchos::RCP<Epetra_CrsMatrix> object, so i think deleting it would solve the problem of memory leak. Any other ideas? I thought maybe the mutliplication result can be stored in the matrix A so that I don't need to allocate a new matrix for result and free the matrix A.

I have to say sorry if the question is very easy but I just couldn't figure it out.

Best regards,
Mehmet Salih YILDIRIM
Hayri Uğur KOLTUK

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


More information about the Trilinos-Users mailing list