[Trilinos-Users] Adding sparse matrices

Williams, Alan B william at sandia.gov
Mon Jun 11 14:52:26 MDT 2007


Maxsim,

Ok, I see what you mean. The error is occurring because Add is
attempting to put a new position into a matrix that has already been
Filled().

Let me look into this issue and I'll get back to you.

Alan


> -----Original Message-----
> From: Maxsim Gibiansky [mailto:maxsim.gibiansky at nist.gov] 
> Sent: Monday, June 11, 2007 2:28 PM
> To: Williams, Alan B
> Cc: trilinos-users at software.sandia.gov
> Subject: RE: [Trilinos-Users] Adding sparse matrices
> 
> > As the function's comments state, it doesn't matter whether 
> the result 
> > has already been 'Filled()' or not, on entry to the function.
> 
> I tried it, and it still gave me an error. Here's code that I 
> can run in the python interpreter:
> 
> from PyTrilinos import Epetra
> from PyTrilinos import EpetraExt
> 
> Comm = Epetra.PyComm()
> Map = Epetra.Map(3,0,Comm)
> A=Epetra.FECrsMatrix(Epetra.Copy, Map,1) 
> B=Epetra.FECrsMatrix(Epetra.Copy, Map,1) 
> C=Epetra.FECrsMatrix(Epetra.Copy, Map,1)
> 
> B[0,0]=1
> C[1,1]=2
> 
> B.FillComplete()
> C.FillComplete()
> 
> EpetraExt.Add(B, False, 1, A, 1)
> EpetraExt.Add(C, False, 1, A, 1)
> 
> This gives me the following error after the last line:
> 
> python2.3: 
> ../../../../packages/epetraext/src/EpetraExt_MatrixMatrix.cpp:1359:
> static int EpetraExt::MatrixMatrix::Add(const 
> Epetra_CrsMatrix&, bool, double, Epetra_CrsMtrix&, double): 
> Assertion `err == 0' failed.
> Aborted
> 
> I'm using trilinos 7.0.8.
> 
> In the description of the Add function, it also says about 
> the FillComplete condition - "If it has, then B's graph must 
> already contain all nonzero locations that will be produced 
> when forming the sum." I think I can only use this if I am 
> sure that the two matrices can be added without generating 
> any new nonzero locations in B, which might not be the case. 
> 
> I guess in this test case, it's possible to get around the 
> problem by just accumulating them in the right order, so that 
> two matrices that are both
> Filled() never need to be added. However, that doesn't seem 
> like it'll work in the general case, since some of the 
> matrices to be added together might already start out 
> FillComplete()d if they were generated by, say, a multiplication. 
> 
> -Maxsim Gibiansky
> 
> Quoting "Williams, Alan B" <william at sandia.gov>:
> 
> > Maxsim,
> > 
> > The function MatrixMatrix::Add(A, transA, scalA, B, scalB) 
> Forms the 
> > result B = scalA * A + scalB * B
> > 
> > As the function's comments state, it doesn't matter whether 
> the result 
> > has already been 'Filled()' or not, on entry to the function.
> > 
> > In your scenario, A is the result matrix which makes the above 
> > somewhat confusing, but the bottom line is, you shouldn't 
> get an error 
> > if the result matrix is already 'Filled()' on entry to the 
> second Add.
> > 
> > Alan
> > 
> > 
> > > -----Original Message-----
> > > From: trilinos-users-bounces at software.sandia.gov
> > > [mailto:trilinos-users-bounces at software.sandia.gov] On Behalf Of 
> > > Maxsim Gibiansky
> > > Sent: Monday, June 11, 2007 1:02 PM
> > > To: trilinos-users at software.sandia.gov
> > > Subject: [Trilinos-Users] Adding sparse matrices
> > > 
> > > I've been using Trilinos (via PyTrilinos), and I have some other 
> > > questions.
> > > 
> > > I'm in a situation where I have sparse matrices
> > > (Epetra.FECrsMatrix) for several terms in an equation; 
> I'd like to 
> > > be able to add them. I've tried using the Add function in 
> > > EpetraExt::MatrixMatrix, but I've had problems - at the 
> end of the 
> > > operation, the destination matrix has
> > > FillComplete() called on it. So if I wanted to use it to do 
> > > something like
> > > 
> > > Add(B, False, 1, A, 1) // A += B
> > > Add(C, False, 1, A, 1) // A += C
> > > 
> > > I might get an error, because A will already be Filled() 
> after the 
> > > first call.
> > > 
> > > Am I missing something? I can't seem to figure out how to add a 
> > > group of matrices together, and it seems like there 
> should be some 
> > > way of doing this.
> > > 
> > > -Maxsim Gibiansky
> > > 
> > > 
> > > _______________________________________________
> > > Trilinos-Users mailing list
> > > Trilinos-Users at software.sandia.gov
> > > http://software.sandia.gov/mailman/listinfo/trilinos-users
> > > 
> > 
> > 
> 
> 
> 
> 
> 
> 



More information about the Trilinos-Users mailing list