[Trilinos-Users] [EXTERNAL] Re: Problem inserting non-local data in Epetra_FECrsMatrix

Heroux, Michael A maherou at sandia.gov
Wed Feb 15 11:48:04 MST 2012


Erhan,

Thanks for this post.  I tried your suggestion and it seems to work using Eric's original code with the following change as you suggest:

      if (A.LRID(globalRowIdx)==-1)
        ierr = A.SumIntoGlobalValues(globalRowIdx, numEntries,
                                    &values[0], &indices[0]);
      else
        ierr = A.InsertGlobalValues(globalRowIdx, numEntries,
                                    &values[0], &indices[0]);

I will talk with Alan Williams (author of this class) about this behavior.  It seems non-intuitive to me that InsertGlobalValues does not work, and does not report any error.

Mike


From: erhan turan <turane at gmail.com<mailto:turane at gmail.com>>
Date: Wed, 15 Feb 2012 19:28:35 +0100
To: Eric Marttila <eric.marttila at thermoanalytics.com<mailto:eric.marttila at thermoanalytics.com>>
Cc: <trilinos-users at software.sandia.gov<mailto:trilinos-users at software.sandia.gov>>
Subject: [EXTERNAL] Re: [Trilinos-Users] Problem inserting non-local data in Epetra_FECrsMatrix

Hello,

I had a similar problem before. My colleague Cyril Flaig found a workaround. For non-local entries, you can try SumIntoGlobalValues rather than InsertGlobalValues. I mean, If A.LRID()== -1 , use SumIntoGlobalValues with global row index and use InsertGlobalValues otherwise. That might work.

For that purpose, what is defined was two pointers say AF and A using Epetra_FECrsMatrix and Epetra_CrsMatrix, respectively. Here, we first declared AF using the map and later just saying A=AF. Later all methods are used for A: A->InsertGlobalValues() and A->SumIntoGlobalValues() . Fill is completed using AF->GlobalAssemble(). If you try, A->FillComplete() it wont fill non-local entries. That might be not the best solution. Maybe there is a simpler way. But it will probably work.

Since you mentioned this issue, I want to make a reference to Epetra_FECrsMatrix documentation. There it reads:

InsertGlobalValues()<http://trilinos.sandia.gov/packages/docs/r10.8/packages/epetra/doc/html/classEpetra__FECrsMatrix.html#a54820e91fb84c7c594c4e6f3ea290fe4> inserts values into the matrix only if the graph has not yet been finalized (FillComplete()<http://trilinos.sandia.gov/packages/docs/r10.8/packages/epetra/doc/html/classEpetra__CrsMatrix.html#a355b4e8d1e28b16bacf37ab637ac0c3f> has not yet been called). For non-local values, the call to InsertGlobalValues()<http://trilinos.sandia.gov/packages/docs/r10.8/packages/epetra/doc/html/classEpetra__FECrsMatrix.html#a54820e91fb84c7c594c4e6f3ea290fe4> may succeed but theGlobalAssemble()<http://trilinos.sandia.gov/packages/docs/r10.8/packages/epetra/doc/html/classEpetra__FECrsMatrix.html#aca5191dcdfb99eb6316b86fee459b243> method may then fail because the non-local data is not actually inserted in the underlying matrix until GlobalAssemble()<http://trilinos.sandia.gov/packages/docs/r10.8/packages/epetra/doc/html/classEpetra__FECrsMatrix.html#aca5191dcdfb99eb6316b86fee459b243> is called.

This might be the reason that FillComplete() alone does not work but those overloaded SumInto methods over GlobalAssemble() work. However I still believe that this sentence may be vague :)

Erhan Turan
Department of Computer Science
ETH Zurich


On Wed, Feb 15, 2012 at 4:19 PM, Eric Marttila <eric.marttila at thermoanalytics.com<mailto:eric.marttila at thermoanalytics.com>> wrote:
Dear all,

I've been using Epetra_CrsMatrix to store a linear system, but recently
encountered a situation where I need to fill the entire matrix from the master
process.  From reading the Epetra documentation I believe I need to switch to
a EpetraFECrs_Matrix to do this.

I changed my code to use the EpetraFECrs_Matrix, but my call to
InsertGlobalValues() fails on the first row of the matrix that is non-local to
the master process.

I've attached a small section of code that reproduces the behavior.  It runs
fine with 1 process, but fails when inserting row index 2 when run with 2
processors.

I would appreciate any thoughts on what I'm missing with my usage of
EpetraFECrs_Matrix.

Thanks.
--Eric

p.s. I'm running Trilinos 10.8.5 on linux x86_64, using mpich2.

--
Eric A. Marttila
ThermoAnalytics, Inc.
23440 Airpark Blvd.
Calumet, MI 49913

email: Eric.Marttila at ThermoAnalytics.com<mailto:Eric.Marttila at ThermoAnalytics.com>
phone: 810-636-2443
fax:   906-482-9755
web: http://www.thermoanalytics.com

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


_______________________________________________ Trilinos-Users mailing list Trilinos-Users at software.sandia.gov<mailto: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/20120215/6bde8c43/attachment-0001.html 


More information about the Trilinos-Users mailing list