[Trilinos-Users] Tpetra alternative to Epetra_MpiComm::SumAll()

Hoemmen, Mark mhoemme at sandia.gov
Tue Apr 21 15:01:36 EDT 2015



On 4/21/15, 12:00 PM, "trilinos-users-request at software.sandia.gov"
<trilinos-users-request at software.sandia.gov> wrote:
>Today's Topics:
>
>   1. Tpetra alternative to Epetra_MpiComm::SumAll() (Truman Ellis)
>
>----------------------------------------------------------------------
>
>Message: 1
>Date: Mon, 20 Apr 2015 20:03:51 +0000
>From: Truman Ellis <truman at ices.utexas.edu>
>To: trilinos-users at software.sandia.gov
>Subject: [Trilinos-Users] Tpetra alternative to
>	Epetra_MpiComm::SumAll()
>Message-ID:
>	<CAPiRBU===cZvugF15aCR-oQFjYo=2Vfk6Q-5iOZVe6sOTyLcEw at mail.gmail.com>
>Content-Type: text/plain; charset="utf-8"
>
>I'm currently in the progress of transitioning our code from explicit use
>of doubles and Epetra to a templated version with Tpetra with the
>intention
>of supporting std::complex<double> solves. In part of our code we have
>Epetra_MpiComm Comm(MPI_COMM_WORLD);
>Comm.SumAll(&mySumCopy, &mySum, 1);
>where mySumCopy and mySum are doubles.
>I noticed there was documentation for Tpetra:Comm with SumAll (I also need
>ScanSum) which matched this Epetra version:
>http://trilinos.org/docs/r8.0/packages/tpetra/doc/html/classTpetra_1_1Comm
>.html
>But when I tried to use this, I got
>Tpetra_MpiComm.hpp: No such file or directory.
>So I'm guessing this was deprecated in favor of all MPI stuff being done
>with Teuchos::Comm, but Teuchos::Comm doesn't have SumAll or ScanSum
>methods. I'm guessing the equivalent calls would be to
>Teuchos::Comm::reduceAllAndScatter and scan with
>Teuchos::SumValueReductionOp passed in as the first argument. But I'm not
>really sure what arguments I need to complete the call to get the same
>results as the Epetra_Comm::SumAll and ScanSum methods produce.

Hi Truman -- Please look in teuchos/core/src/Teuchos_CommHelpers.hpp for
various overloads of the Teuchos::reduceAll function.  For example:

using Teuchos::Comm;
using Teuchos::outArg;

using Teuchos::RCP;
using Teuchos::REDUCE_SUM;

using Teuchos::reduceAll;

RCP<const Comm<int> > comm = ...;

double localSum = ...;
double globalSum = 0.0;

// Sum up the localSum values over all processes in comm,
// and write the result to globalSum on all processes.
reduceAll (*comm, REDUCE_SUM, localSum, outArg (globalSum));

You probably want the variants of reduceAll that take an EReductionType
enum value.  Other values are REDUCE_MIN and REDUCE_MAX, which do what you
think.  If you have an array of values to reduce, you can use the
five-argument version of reduceAll:

double* localSums = new double [count];
// ... compute local sums ...
double* globalSums = new double [count];
reduceAll (*comm, REDUCE_SUM, count, localSums, globalSums);

Hope this helps!
mfh



More information about the Trilinos-Users mailing list