46 #ifndef MUELU_AMGXOPERATOR_DECL_HPP 
   47 #define MUELU_AMGXOPERATOR_DECL_HPP 
   49 #if defined (HAVE_MUELU_AMGX) 
   52 #include <Tpetra_Operator.hpp> 
   53 #include <Tpetra_CrsMatrix.hpp> 
   54 #include <Tpetra_MultiVector.hpp> 
   55 #include <Tpetra_Distributor.hpp> 
   56 #include <Tpetra_HashTable.hpp> 
   57 #include <Tpetra_Import.hpp> 
   58 #include <Tpetra_Import_Util.hpp> 
   62 #include "MueLu_TpetraOperator.hpp" 
   65 #include <cuda_runtime.h> 
   81   class AMGXOperator : 
public TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>, 
public BaseClass {
 
   88     typedef Tpetra::Map<LO,GO,NO>            
Map;
 
  150     typedef Tpetra::Map<LO,GO,NO>            
Map;
 
  155                    const int* nbrs, 
const Map& map, 
const std::string& label) {
 
  156       for (
int p = 0; p < comm->getSize(); p++) {
 
  157         if (comm->getRank() == p) {
 
  158           std::cout << 
"========\n" << label << 
", lid (gid), PID  " << p << 
"\n========" << std::endl;
 
  160           for (
size_t i = 0; i < vec.size(); ++i) {
 
  161             std::cout << 
"   neighbor " << nbrs[i] << 
" :";
 
  162             for (
size_t j = 0; j < vec[i].size(); ++j)
 
  163               std::cout << 
" " << vec[i][j] << 
" (" << map.getGlobalElement(perm[vec[i][j]]) << 
")";
 
  164             std::cout << std::endl;
 
  166           std::cout << std::endl;
 
  180       int numProcs = comm->getSize();
 
  181       int myRank   = comm->getRank();
 
  193       AMGX_SAFE_CALL(AMGX_install_signal_handler());
 
  196         AMGX_SAFE_CALL(AMGX_config_create_from_file(&Config_, (
const char *) &configs.
get<std::string>(
"json file")[0]));
 
  198         std::ostringstream oss;
 
  201         for (itr = configs.
begin(); itr != configs.
end(); ++itr) {
 
  202           const std::string&    name  = configs.
name(itr);
 
  204           oss << name << 
"=" << filterValueToString(entry) << 
", ";
 
  207         std::string configString = oss.str();
 
  208         if (configString == 
"") {
 
  212         AMGX_SAFE_CALL(AMGX_config_create(&Config_, configString.c_str()));
 
  228       MPI_Comm mpiComm = *rawMpiComm;
 
  233         AMGX_resources_create_simple(&Resources_, Config_);
 
  237         cudaGetDeviceCount(&numGPUDevices);
 
  238         int device[] = {(comm->getRank() % numGPUDevices)};
 
  240         AMGX_config_add_parameters(&Config_, 
"communicator=MPI");
 
  242         AMGX_resources_create(&Resources_, Config_, &mpiComm, 1, device);
 
  244         AMGX_resources_create(&Resources_, Config_, MPI_COMM_WORLD, 1, device);
 
  248       AMGX_Mode mode = AMGX_mode_dDDI;
 
  249       AMGX_solver_create(&Solver_, Resources_, mode,  Config_);
 
  250       AMGX_matrix_create(&A_,      Resources_, mode);
 
  251       AMGX_vector_create(&X_,      Resources_, mode);
 
  252       AMGX_vector_create(&Y_,      Resources_, mode);
 
  257       std::vector<int> amgx2muelu;
 
  265         Tpetra::Distributor distributor = importer->getDistributor();
 
  267         Array<int> sendRanks = distributor.getProcsTo();
 
  268         Array<int> recvRanks = distributor.getProcsFrom();
 
  270         std::sort(sendRanks.
begin(), sendRanks.
end());
 
  271         std::sort(recvRanks.
begin(), recvRanks.
end());
 
  274         if (sendRanks.
size() != recvRanks.
size()) {
 
  277           for (
int i = 0; i < sendRanks.
size(); i++) {
 
  278             if (recvRanks[i] != sendRanks[i])
 
  284                                    "This is not the case: we send to {" << sendRanks << 
"} and receive from {" << recvRanks << 
"}");
 
  286         int        num_neighbors = sendRanks.
size();  
 
  287         const int* neighbors     = &sendRanks[0];
 
  293         Tpetra::Details::HashTable<int,int> hashTable(3*num_neighbors);
 
  294         for (
int i = 0; i < num_neighbors; i++)
 
  295           hashTable.add(neighbors[i], i);
 
  301         Tpetra::Import_Util::getPids(*importer, importPIDs, 
true);
 
  307         int N = rowMap->getNodeNumElements(), Nc = colMap->getNodeNumElements();
 
  308         muelu2amgx_.resize(Nc, -1);
 
  310         int numUniqExports = 0;
 
  311         for (
int i = 0; i < exportLIDs.
size(); i++)
 
  312           if (muelu2amgx_[exportLIDs[i]] == -1) {
 
  314             muelu2amgx_[exportLIDs[i]] = -2;
 
  317         int localOffset = 0, exportOffset = N - numUniqExports;
 
  319         for (
int i = 0; i < exportLIDs.
size(); i++)
 
  320           if (muelu2amgx_[exportLIDs[i]] < 0) 
 
  321             muelu2amgx_[exportLIDs[i]] = exportOffset++;
 
  323         for (
int i = 0; i < N; i++)
 
  324           if (muelu2amgx_[i] == -1)
 
  325             muelu2amgx_[i] = localOffset++;
 
  327         int importOffset = N;
 
  328         for (
int k = 0; k < num_neighbors; k++)
 
  329           for (
int i = 0; i < importPIDs.
size(); i++)
 
  330             if (importPIDs[i] != -1 && hashTable.get(importPIDs[i]) == k)
 
  331               muelu2amgx_[i] = importOffset++;
 
  333         amgx2muelu.resize(muelu2amgx_.size());
 
  334         for (
int i = 0; i < (int)muelu2amgx_.size(); i++)
 
  335           amgx2muelu[muelu2amgx_[i]] = i;
 
  338         std::vector<std::vector<int> > sendDatas (num_neighbors);
 
  339         std::vector<int>               send_sizes(num_neighbors, 0);
 
  340         for (
int i = 0; i < exportPIDs.
size(); i++) {
 
  341           int index = hashTable.get(exportPIDs[i]);
 
  342           sendDatas [index].push_back(muelu2amgx_[exportLIDs[i]]);
 
  347         std::vector<const int*> send_maps(num_neighbors);
 
  348         for (
int i = 0; i < num_neighbors; i++)
 
  349           send_maps[i] = &(sendDatas[i][0]);
 
  355         std::vector<std::vector<int> > recvDatas (num_neighbors);
 
  356         std::vector<int>               recv_sizes(num_neighbors, 0);
 
  357         for (
int i = 0; i < importPIDs.
size(); i++)
 
  358           if (importPIDs[i] != -1) {
 
  359             int index = hashTable.get(importPIDs[i]);
 
  360             recvDatas [index].push_back(muelu2amgx_[i]);
 
  365         std::vector<const int*> recv_maps(num_neighbors);
 
  366         for (
int i = 0; i < num_neighbors; i++)
 
  367           recv_maps[i] = &(recvDatas[i][0]);
 
  372         AMGX_SAFE_CALL(AMGX_matrix_comm_from_maps_one_ring(A_, 1, num_neighbors, neighbors, &send_sizes[0], &send_maps[0], &recv_sizes[0], &recv_maps[0]));
 
  374         AMGX_vector_bind(X_, A_);
 
  375         AMGX_vector_bind(Y_, A_);
 
  379       matrixTransformTimer->
start();
 
  384       inA->getAllValues(ia_s, ja, a);
 
  387       for (
int i = 0; i < ia.size(); i++)
 
  388         ia[i] = Teuchos::as<int>(ia_s[i]);
 
  390       N_      = inA->getNodeNumRows();
 
  391       int nnz = inA->getNodeNumEntries();
 
  393       matrixTransformTimer->
stop();
 
  400       matrixTimer->
start();
 
  402         AMGX_matrix_upload_all(A_, N_, nnz, 1, 1, &ia[0], &ja[0], &a[0], NULL);
 
  406         std::vector<int>    ia_new(ia.size());
 
  407         std::vector<int>    ja_new(ja.
size());
 
  408         std::vector<double> a_new (a.
size());
 
  411         for (
int i = 0; i < N_; i++) {
 
  412           int oldRow = amgx2muelu[i];
 
  414           ia_new[i+1] = ia_new[i] + (ia[oldRow+1] - ia[oldRow]);
 
  416           for (
int j = ia[oldRow]; j < ia[oldRow+1]; j++) {
 
  417             int offset = j - ia[oldRow];
 
  418             ja_new[ia_new[i] + offset] = muelu2amgx_[ja[j]];
 
  419             a_new [ia_new[i] + offset] = a[j];
 
  427             for (
int j = ia_new[i]; j < ia_new[i+1]-1; j++)
 
  428               if (ja_new[j] > ja_new[j+1]) {
 
  429                 std::swap(ja_new[j], ja_new[j+1]);
 
  430                 std::swap(a_new [j], a_new [j+1]);
 
  433           } 
while (swapped == 
true);
 
  436         AMGX_matrix_upload_all(A_, N_, nnz, 1, 1, &ia_new[0], &ja_new[0], &a_new[0], NULL);
 
  441       domainMap_ = inA->getDomainMap();
 
  442       rangeMap_  = inA->getRangeMap();
 
  445       realSetupTimer->
start();
 
  446       AMGX_solver_setup(Solver_, A_);
 
  447       realSetupTimer->
stop();
 
  459       AMGX_SAFE_CALL(AMGX_solver_destroy(Solver_));
 
  460       AMGX_SAFE_CALL(AMGX_vector_destroy(X_));
 
  461       AMGX_SAFE_CALL(AMGX_vector_destroy(Y_));
 
  462       AMGX_SAFE_CALL(AMGX_matrix_destroy(A_));
 
  463       AMGX_SAFE_CALL(AMGX_resources_destroy(Resources_));
 
  464       AMGX_SAFE_CALL(AMGX_config_destroy(Config_));
 
  495       AMGX_matrix_get_size(A_, &n, &sizeX, &sizeY);
 
  501       AMGX_solver_get_iterations_number(Solver_, &it);
 
  506       AMGX_SOLVE_STATUS status;
 
  507       AMGX_solver_get_status(Solver_, &status);
 
  516     AMGX_matrix_handle      
A_;
 
  517     AMGX_vector_handle      
X_;
 
  518     AMGX_vector_handle      
Y_;
 
  533 #endif //HAVE_MUELU_AMGX 
  534 #endif // MUELU_AMGXOPERATOR_DECL_HPP 
RCP< MueLu::Hierarchy< SC, LO, GO, NO > > GetHierarchy() const 
 
RCP< const Map > domainMap_
 
const std::string & name() const 
 
virtual ~AMGXOperator()
Destructor. 
 
ConstIterator end() const 
 
MueLu::DefaultLocalOrdinal LocalOrdinal
 
std::string toString(const T &what)
Little helper function to convert non-string types to strings. 
 
AMGXOperator(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &inA, Teuchos::ParameterList ¶mListIn)
 
Tpetra::Map< LO, GO, NO > Map
 
T & get(const std::string &name, T def_value)
 
AMGX_config_handle Config_
 
void printMaps(Teuchos::RCP< const Teuchos::Comm< int > > &comm, const std::vector< std::vector< int > > &vec, const std::vector< int > &perm, const int *nbrs, const Map &map, const std::string &label)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
Teuchos::RCP< const Map > getRangeMap() const 
Returns the Tpetra::Map object associated with the range of this operator. 
 
AMGX_solver_handle Solver_
 
std::vector< int > muelu2amgx_
 
virtual ~AMGXOperator()
Destructor. 
 
Tpetra::Map< LO, GO, NO > Map
 
RCP< const Map > rangeMap_
 
RCP< Teuchos::Time > solverTimer_
 
static RCP< Time > getNewTimer(const std::string &name)
 
bool isParameter(const std::string &name) const 
 
void start(bool reset=false)
 
MueLu::DefaultScalar Scalar
 
MueLu::DefaultGlobalOrdinal GlobalOrdinal
 
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const 
Returns a solution for the linear system AX=Y in the Tpetra::MultiVector X. 
 
params_t::ConstIterator ConstIterator
 
ConstIterator begin() const 
 
std::string filterValueToString(const Teuchos::ParameterEntry &entry)
 
const ParameterEntry & entry(ConstIterator i) const 
 
any & getAny(bool activeQry=true)
 
Teuchos::RCP< const Map > getDomainMap() const 
Returns the Tpetra::Map object associated with the domain of this operator. 
 
AMGXOperator(const Teuchos::RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > &InA, Teuchos::ParameterList ¶mListIn)
Constructor. 
 
AMGX_resources_handle Resources_
 
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator. 
 
Tpetra::MultiVector< SC, LO, GO, NO > MultiVector
 
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
 
Exception throws to report errors in the internal logical of the program. 
 
Tpetra::MultiVector< SC, LO, GO, NO > MultiVector
 
RCP< MueLu::Hierarchy< SC, LO, GO, NO > > GetHierarchy() const 
 
RCP< Teuchos::Time > vectorTimer1_
 
bool hasTransposeApply() const 
Indicates whether this operator supports applying the adjoint operator. 
 
RCP< Teuchos::Time > vectorTimer2_
 
AMGX_SOLVE_STATUS getStatus()
 
Adapter for AmgX library from Nvidia.