20 #ifndef CODE_MPILIB_MPINODE_HPP_
21 #define CODE_MPILIB_MPINODE_HPP_
31 template<
class Weight,
class NodeDistribution>
36 const NodeDistribution& nodeDistribution,
38 const std::string& name) :
39 _pAlgorithm(algorithm.clone()),
42 _rLocalNodes(localNode),
43 _rNodeDistribution(nodeDistribution),
47 template<
class Weight,
class NodeDistribution>
51 template<
class Weight,
class NodeDistribution>
63 _pAlgorithm->evolveNodeState(_precursorActivity, _weights, time,
65 Time t_ret = _pAlgorithm->getCurrentTime();
72 this->setActivity(_pAlgorithm->getCurrentRate());
77 return _pAlgorithm->getCurrentTime();
80 template<
class Weight,
class NodeDistribution>
82 _pAlgorithm->prepareEvolve(_precursorActivity, _weights, _precursorTypes);
86 template<
class Weight,
class NodeDistribution>
91 _pAlgorithm->configure(simParam);
94 this->setActivity(_pAlgorithm->getCurrentRate());
96 _pHandler = std::shared_ptr<report::handler::AbstractReportHandler>(
99 _pHandler->initializeHandler(_nodeId);
102 template<
class Weight,
class NodeDistribution>
104 const Weight& weight,
NodeType nodeType) {
105 _precursors.push_back(nodeId);
106 _precursorTypes.push_back(nodeType);
107 _weights.push_back(weight);
109 _precursorActivity.resize(_precursors.size());
112 template<
class Weight,
class NodeDistribution>
114 _successors.push_back(nodeId);
117 template<
class Weight,
class NodeDistribution>
122 template<
class Weight,
class NodeDistribution>
124 _activity = activity;
127 template<
class Weight,
class NodeDistribution>
133 template<
class Weight,
class NodeDistribution>
137 for (
auto it = _precursors.begin(); it != _precursors.end(); it++, i++) {
139 if (_rNodeDistribution.isLocalNode(*it)) {
140 _precursorActivity[
i] =
141 _rLocalNodes.find(*it)->second.getActivity();
145 _precursorActivity[
i]);
150 template<
class Weight,
class NodeDistribution>
153 for (
auto& it : _successors) {
155 if (!_rNodeDistribution.isLocalNode(it)) {
162 template<
class Weight,
class NodeDistribution>
165 std::vector<report::ReportValue> vec_values;
169 Rate(this->getActivity()), this->_nodeId,
170 _pAlgorithm->getGrid(this->_nodeId,
false), type, vec_values, _rLocalNodes.size());
171 _pHandler->writeReport(report);
177 Rate(this->getActivity()), this->_nodeId,
178 _pAlgorithm->getGrid(this->_nodeId,_pHandler->isStateWriteMandatory()), type, vec_values, _rLocalNodes.size());
179 _pHandler->writeReport(report);
183 template<
class Weight,
class NodeDistribution>
186 _pHandler->detachHandler(_nodeId);
189 template<
class Weight,
class NodeDistribution>
NodeType
This determines an MPINode's type, which will be checked when Dale's law is set.
Class for nodes in an MPINetwork.
ActivityType getActivity() const
void irecv(int, int, T &) const
A Report is sent by a MPINode when it is queried.
MPINode(const AlgorithmInterface< Weight > &algorithm, NodeType nodeType, NodeId nodeId, const NodeDistribution &nodeDistribution, const std::map< NodeId, MPINode< Weight, NodeDistribution >> &localNode, const std::string &name="")
NodeType getNodeType() const
void configureSimulationRun(const SimulationRunParameter &simParam)
Number getMaximumNumberIterations() const
const double ALGORITHM_NETWORK_DISCREPANCY
Individual algorithms may make small errors in their time keeping compared to the network...
const report::handler::AbstractReportHandler & getHandler() const
virtual AbstractReportHandler * clone() const =0
void addSuccessor(NodeId nodeId)
Parameter determining how a simulation is run. Specifiying begin and end time, log file names...
void reportAll(report::ReportType type) const
void setActivity(ActivityType activity)
void addPrecursor(NodeId nodeId, const Weight &weight, NodeType nodeType)
The objective is find a numerical solution for this equation This requires a numerical representation of the density We will work in the state space of a two dimensional and define a mesh there We first give two examples and then define the general procedure and given in Table for given fixed Delta g see Fig and we will denote coordinates in this dimension by a small letter $v The second dimension can be used to represent parameters as varied as and will represented by $w A strip is constructed by choosing two neighbouring points in state e g and integrating the vector field for a time $T that is assumed to be an integer multiple of a period of time Delta which we assume to be a defining characteristic of the grid Let then the set of points the set of points which is quadrilateral in shape The quadrilateral should be but not necessarily as long as they are but it is convenient to number them in order of creation In the we will assume that strip numbers created by the integration procedure start and are so that the numbers i in each identify a unique strip Strip no is reserved for stationary points There may or more cells in strip The number of cells in strip $i denoted by with $i the strip number and $j the cell as the i
void isend(int, int, const T &) const