[Trilinos-Users] Tpetra::Vector: elementwise operations

Baker, Christopher G. bakercg at ornl.gov
Mon Jan 25 09:21:30 MST 2010


Actually, the get1dView() method returns the local values, so that the indexing of inVectorView[] by global IDs is not correct. get1dView() is, in fact, the getLocalValues() method you were looking for. There is no getGlobalValues() method. One would need to explicitly import values into a new Vector object, then use get1dView() on it.

Furthermore, and I keep meaning to add a note about this, the method replaceLocalValue() is not an efficient method anymore. On GPU platforms, the performance will be horrible. On CPU platforms, there will be some overhead to pay. The reason is because each call to replaceLocalValue() will call get1dViewNonConst(); therefore, it is better to call get1dViewNonConst() up front and use a view of the out vector, as you are currently doing with the in vector.

With Tpetra objects, as with Epetra objects, you very seldom have access to global objects without explicitly requesting that they are transported to your process. Unlike Epetra objects, you don't always have efficient access to local values, because of the potential that the values are stored "far away" (e.g., on a GPU device). Some methods from Epetra_Vector have been moved over, unfortunately without noting possible inefficiency (e.g., replaceLocalValue()). Other methods have been neglected, because of possible inefficiency (e.g., operator[]). One outstanding task for the Tpetra developers is to decide whether we want to support these inefficient methods (I think we do), then to add them, along with a warning that they are not efficient, and a pointer on which approach is better. Some of these methods are very convenient for debugging, and I think their utility there and for novice users (per our emails last week) justifies their existence/maintenance.

Underlying the MultiVector and Vector objects is a mechanism for generically applying elementwise transformations and reductions. We are currently exploring how to expose this for Tpetra users. Such an interface would allow you to write custom reduction/transformation ops (RTOps, a la Trilinos/Thyra) and execute them for/on a Vector. Therefore, users could write custom operations and even have them executed on the GPU in GPU-based platforms, without having to pull data over to the host view methods like get1dView(). I would appreciate your comments on this.

Chris


On 1/24/10 8:51 AM, "Nico Schlömer" <nico.schloemer at ua.ac.be> wrote:

Hi,

thinking about Tpetra vectors, I was asking myself if there's any better
way to do elementwise operations (e.g., sqrt(), conj(),...) than what I do
at the moment:

===================== *snip* =====================
Teuchos::ArrayRCP<const double> inVectorView = inVector.get1dView();
for (int k=0; k<NumLocalElements; k++ ) {
   int globalIndex = Map->getGlobalElement(k);
   outVector->replaceLocalValue( k, 0, inVectorView[globalIndex] );
}
===================== *snap* =====================

As far as I can see, there's no such method as getLocalValue{s}(). Am I
right?

Cheers,
Nico


_______________________________________________
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