Line data Source code
1 : // Copyright (C) 2015 Technische Universitaet Muenchen 2 : // This file is part of the Mamico project. For conditions of distribution 3 : // and use, please see the copyright notice in Mamico's main folder, or at 4 : // www5.in.tum.de/mamico 5 : #ifndef _MOLECULARDYNAMICS_COUPLING_SENDRECV_DATAEXCHANGEFROMMACRO2MD_H_ 6 : #define _MOLECULARDYNAMICS_COUPLING_SENDRECV_DATAEXCHANGEFROMMACRO2MD_H_ 7 : 8 : #include "coupling/CouplingMDDefinitions.h" 9 : #include "coupling/datastructures/CouplingCell.h" 10 : #include "coupling/indexing/IndexingService.h" 11 : #include "coupling/interface/MacroscopicSolverInterface.h" 12 : #include "coupling/sendrecv/DataExchange.h" 13 : 14 : namespace coupling { 15 : namespace sendrecv { 16 : template <unsigned int dim> class DataExchangeFromMacro2MD; 17 : } 18 : } // namespace coupling 19 : 20 : /** data exchange from the macroscopic solver to the MD solver, that is to the 21 : *coupling tool. We transfer the buffers microscopicMass and microscopicMomentum 22 : *of the CouplingCell. The target ranks are determined by the 23 : *getRanksForCouplingCell(): since coupling 24 : *cells may exist on different ranks (due to ghost layers), this method 25 : *determines all ranks with a particular coupling cell and returns a vector 26 : *with all required ranks. The 27 : *source ranks are determined via the macroscopic solver interface's method 28 : *getRanks() 29 : * @brief data exchange from the macroscopic solver to the MD solver. 30 : *Derived from the class coupling::sendrecv::DataExchange 31 : * @tparam dim Number of dimensions; it can be 1, 2 or 3 32 : * @sa DataExchangeFromMD2Macro 33 : * @author Philipp Neumann 34 : */ 35 : template <unsigned int dim> 36 : class coupling::sendrecv::DataExchangeFromMacro2MD : public coupling::sendrecv::DataExchange<coupling::datastructures::CouplingCell<dim>, dim> { 37 : 38 : public: 39 : /** Constructor 40 : * @param interface macroscopic solver interface 41 : * @param tagoffset 0 per default 42 : */ 43 4 : DataExchangeFromMacro2MD(coupling::interface::MacroscopicSolverInterface<dim>* interface, unsigned int topologyOffset, unsigned int tagoffset = 0) 44 4 : : coupling::sendrecv::DataExchange<coupling::datastructures::CouplingCell<dim>, dim>(TAG_FROM_MACRO2MD + tagoffset), _msi(interface), 45 4 : _topologyOffset(topologyOffset) { 46 : #if (COUPLING_MD_DEBUG == COUPLING_MD_YES) 47 : std::cout << "DataExchangeFromMacro2MD initialised..." << std::endl; 48 : #endif 49 : } 50 : /** Destructor */ 51 0 : virtual ~DataExchangeFromMacro2MD() {} 52 : 53 : /** returns the ranks to which a particular cell (at index idx) 54 : *should be sent. 55 : * @param idx 56 : * @return the corresponding ranks, if we need 57 : *information on MD side, otherwise empty vector 58 : */ 59 0 : std::vector<unsigned int> getTargetRanks(I01 idx) override { 60 : // if we need information on MD side, return the respective ranks via 61 0 : if (I08::contains(idx) && !I12::contains(idx)) { 62 0 : return IDXS.getRanksForGlobalIndex(idx, _topologyOffset); 63 : // otherwise return empty vector 64 : } 65 0 : return std::vector<unsigned int>(); 66 : } 67 : 68 : /** returns all ranks from which a particular cell (at index idx) 69 : *is sent. 70 : * @param idx 71 : * @return the corresponding ranks via MacroscopicSolverInterface, if we 72 : *need information on MD side, otherwise empty vector 73 : */ 74 0 : std::vector<unsigned int> getSourceRanks(I01 idx) override { 75 0 : if (I08::contains(idx) && !I12::contains(idx)) { // Global, no ghost, no MD2Macro 76 0 : return _msi->getSourceRanks(idx); 77 : } 78 0 : return std::vector<unsigned int>(); 79 : } 80 : 81 : /** local rule to read from coupling cell and write data to (e.g. send) 82 : * buffer. We only send the macroscopic mass and macroscopic momentum from MD 83 : * to the macroscopic solver. 84 : * @param buffer 85 : * @param cell 86 : */ 87 0 : void readFromCell(double* const buffer, const coupling::datastructures::CouplingCell<dim>& cell) override { 88 0 : buffer[0] = cell.getMicroscopicMass(); 89 0 : for (unsigned int d = 0; d < dim; d++) { 90 0 : buffer[d + 1] = cell.getMicroscopicMomentum()[d]; 91 : } 92 0 : } 93 : 94 : /** local rule to read from receive buffer and write data to coupling cell 95 : * @param buffer 96 : * @param cell 97 : */ 98 0 : void writeToCell(const double* const buffer, coupling::datastructures::CouplingCell<dim>& cell) override { 99 0 : tarch::la::Vector<dim, double> microscopicMomentum(0.0); 100 0 : for (unsigned int d = 0; d < dim; d++) { 101 0 : microscopicMomentum[d] = buffer[1 + d]; 102 : } 103 0 : cell.setMicroscopicMomentum(microscopicMomentum); 104 0 : cell.setMicroscopicMass(buffer[0]); 105 0 : } 106 : 107 : /** returns the number of doubles that are sent per coupling cell. @return 108 : * 1+dim */ 109 0 : unsigned int getDoublesPerCell() const override { 110 : // 1 double: microscopic mass; dim doubles: microscopic momentum 111 0 : return 1 + dim; 112 : } 113 : 114 : private: 115 : coupling::interface::MacroscopicSolverInterface<dim>* _msi; 116 : unsigned int _topologyOffset; 117 : }; 118 : #endif // _MOLECULARDYNAMICS_COUPLING_SENDRECV_DATAEXCHANGEFROMMACRO2MD_H_