CoDiPack  2.2.0
A Code Differentiation Package
SciComp TU Kaiserslautern
Loading...
Searching...
No Matches
Example 13 - MPI communication

Goal: Learn how MPI communication is differentiated with CoDiPack and MeDiPack.

Prerequisite: Tutorial 2 - Reverse mode AD

Full code:

#include <codi.hpp>
#include <medi/medi.hpp>
#include <iostream>
using namespace medi;
#include "codi/tools/mpi/codiMpiTypes.hpp"
using Tape = typename Real::Tape;
using MpiTypes = codi::CoDiMpiTypes<Real>;
MpiTypes* mpiTypes;
int main(int nargs, char** args) {
AMPI_Init(&nargs, &args); // Step 1: Replace all MPI_* functions and types with AMPI_*
mpiTypes = new MpiTypes(); // Step 2: Create the CoDiPack MPI types.
int rank;
int size;
AMPI_Comm_rank(AMPI_COMM_WORLD, &rank);
AMPI_Comm_size(AMPI_COMM_WORLD, &size);
if(size != 2) {
std::cout << "Please start the tutorial with two processes." << std::endl;
} else {
Tape& tape = Real::getTape();
tape.setActive();
Real a = 3.0;
if( 0 == rank ) {
tape.registerInput(a);
AMPI_Send(&a, 1, mpiTypes->MPI_TYPE, 1, 42, AMPI_COMM_WORLD); // Step 3: Use the CoDiPack MPI type as the data type.
} else {
AMPI_Recv(&a, 1, mpiTypes->MPI_TYPE, 0, 42, AMPI_COMM_WORLD, AMPI_STATUS_IGNORE);
tape.registerOutput(a);
a.setGradient(100.0);
}
tape.setPassive();
tape.evaluate();
if(0 == rank) {
std::cout << "Adjoint of 'a' on rank 0 is: " << a.getGradient() << std::endl;
}
}
delete mpiTypes; // Step 4: Cleanup the created CoDiPack MPI types.
AMPI_Finalize();
}
#include <medi/medi.cpp>
RealReverseGen< double > RealReverse
Definition: codi.hpp:120
Represents a concrete lvalue in the CoDiPack expression tree.
Definition: activeType.hpp:52
static Tape & getTape()
Get a reference to the tape which manages this expression.
Definition: activeType.hpp:99
T_Tape Tape
See ActiveType.
Definition: activeType.hpp:55
MPI datatype implementation for CoDipack types in the type wrapper of MeDiPack.
Definition: codiMpiTypes.hpp:62
void setGradient(Gradient const &g)
Set the gradient of this lvalue in the tape.
Definition: lhsExpressionInterface.hpp:120
Gradient getGradient() const
Get the gradient of this lvalue from the tape.
Definition: lhsExpressionInterface.hpp:115

The differentiation of MPI communication is done via the MeDiPack library (https://www.scicomp.uni-kl.de/software/medi). As the example shows all 'MPI_*" routines and functions need to be replaced with the MeDiPack wrapper functions 'AMPI_*'. Otherwise the MPI functions can be used as usual. The only difference is, that the CoDiPack MPI Datatype needs to be used in all communications where CoDiPack types are involved.

If different CoDiPack types are used (e.g. codi::RealReverse and codi::RealReverseIndex), then for each one a CoDiPack MPI Datatype needs to be created.