[Trilinos-Users] Defining bordered problems efficiently: Epetra_Import

Heroux, Michael A maherou at sandia.gov
Tue Nov 25 12:55:42 MST 2008


Daniele,

If I understand what you are doing, then one way to eliminate the extra copy of x values is to use View mode with the Epetra_Vector objects.  You will still have overlappedXExtMap and xExtMap, and the corresponding vectors overlappedXExt and xExt.  Furthermore, you would create overlappedXExt in the standard way.  However, when you create xExt, you should pass in the pointer to the first value in overlappedXExt using View mode.

If you do this, Epetra's import function is smart enough to detect that all but the last value are the same address and it will not do a redundant copy.

A further optimization would be to avoid the import, since you are really doing a broadcast.  You could use the Broadcast() method on Epetra_Comm to send to all processors and pass a pointer to the end of  overlappedXExt vector as the target buffer for the broadcast.

I hope this makes sense.

Mike


On 11/24/08 6:53 AM, "Daniele Avitabile" <d.avitabile at gmail.com> wrote:

Hi everybody,

the solution to the PDEs I am solving is stored in an Epetra_Vector x, based on a certain linear Epetra_Map xMap. Now, I need to add an extra ODE to my system, defined on the boundary. At the present stage, whenever I want to compute the right-hand side of my PDE, I need to access a new "bordered" Epetra_Vector xExt, based on a new map Epetra_Map xExtMap.

The problem is that, from each processor, I need access to the last component of xExt, the one on the boundary.

The way I am doing this is
1. define an Epetra_Map overlappedXExtMap. All the GID of xExtMap are also contained in overlappedXExtMap, plus the last component of xExt is replicated on all processors.
2. define an importer with xExtMap as source and overlappedXExtMap as target
3. each time I evaluate the right-hand side of my equation, use the importer.

I rely on the importer massively, because I also use it for printing purposes: I want to separate the solution in the interior of my PDE from the additional scalar component on the boundary.

My questions:
a. do you have a more efficient way to do this? It seems to me that I am copying a lot of unnecessary data, namely all the components relative to the GID of xEtxMap that are also in overlappedXExtMap.
b. how does the importer Work? Does he really copy all the elements?

Thanks in advance.

Daniele


-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20081125/caaddb10/attachment.html 


More information about the Trilinos-Users mailing list