Program Listing for File MessageArrayDevice.cuh
↰ Return to documentation for file (include/flamegpu/runtime/messaging/MessageArray/MessageArrayDevice.cuh
)
#ifndef INCLUDE_FLAMEGPU_RUNTIME_MESSAGING_MESSAGEARRAY_MESSAGEARRAYDEVICE_CUH_
#define INCLUDE_FLAMEGPU_RUNTIME_MESSAGING_MESSAGEARRAY_MESSAGEARRAYDEVICE_CUH_
#include "flamegpu/runtime/messaging/MessageArray.h"
#include "flamegpu/runtime/messaging/MessageBruteForce/MessageBruteForceDevice.cuh"
namespace flamegpu {
class MessageArray::In {
friend class Message;
public:
class Message {
const MessageArray::In &_parent;
size_type index;
public:
__device__ Message(const MessageArray::In &parent, const size_type _index) : _parent(parent), index(_index) {}
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
__device__ Message(const MessageArray::In& parent) : _parent(parent), index(0) {}
#endif
__device__ bool operator==(const Message& rhs) const { return this->index == rhs.index; }
__device__ bool operator!=(const Message& rhs) const { return this->index != rhs.index; }
__device__ size_type getIndex() const { return this->index; }
template<typename T, unsigned int N>
__device__ T getVariable(const char(&variable_name)[N]) const;
template<typename T, flamegpu::size_type N, unsigned int M>
__device__ T getVariable(const char(&variable_name)[M], unsigned int index) const;
};
class WrapFilter {
friend class Message;
public:
class Message {
const WrapFilter&_parent;
int relative_cell;
size_type index_1d = 0;
public:
__device__ Message(const WrapFilter&parent, const int relative_x)
: _parent(parent) {
relative_cell = relative_x;
}
__device__ bool operator==(const Message& rhs) const {
return this->relative_cell == rhs.relative_cell;
// && this->_parent.loc == rhs._parent.loc;
}
__device__ bool operator!=(const Message& rhs) const { return !(*this == rhs); }
__device__ inline Message& operator++();
__device__ size_type getX() const {
return (this->_parent.loc + relative_cell + this->_parent.length) % this->_parent.length;
}
__device__ size_type getOffsetX() const {
return relative_cell;
}
template<typename T, unsigned int N>
__device__ T getVariable(const char(&variable_name)[N]) const;
template<typename T, flamegpu::size_type N, unsigned int M> __device__
T getVariable(const char(&variable_name)[M], unsigned int index) const;
};
class iterator {
Message _message;
public:
__device__ iterator(const WrapFilter&parent, const int relative_x)
: _message(parent, relative_x) {
// Increment to find first message
++_message;
}
__device__ iterator& operator++() { ++_message; return *this; }
__device__ iterator operator++(int) {
iterator temp = *this;
++*this;
return temp;
}
__device__ bool operator==(const iterator& rhs) const { return _message == rhs._message; }
__device__ bool operator!=(const iterator& rhs) const { return _message != rhs._message; }
__device__ Message& operator*() { return _message; }
__device__ Message* operator->() { return &_message; }
};
inline __device__ WrapFilter(const size_type _length, const size_type x, const size_type _radius);
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
inline __device__ WrapFilter();
#endif
inline __device__ iterator begin(void) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (!this->length)
return iterator(*this, radius);
#endif
// Bin before initial bin, as the constructor calls increment operator
return iterator(*this, -static_cast<int>(radius) - 1);
}
inline __device__ iterator end(void) const {
// Final bin, as the constructor calls increment operator
return iterator(*this, radius);
}
private:
size_type loc;
const size_type radius;
const size_type length;
};
class Filter {
friend class Message;
public:
class Message {
const Filter& _parent;
int relative_cell;
size_type index_1d = 0;
public:
__device__ Message(const Filter& parent, const int relative_x)
: _parent(parent) {
relative_cell = relative_x;
}
__device__ bool operator==(const Message& rhs) const {
return this->relative_cell == rhs.relative_cell;
// && this->_parent.loc == rhs._parent.loc;
}
__device__ bool operator!=(const Message& rhs) const { return !(*this == rhs); }
__device__ inline Message& operator++();
__device__ size_type getX() const {
return this->_parent.loc + relative_cell;
}
__device__ int getOffsetX() const {
return this->relative_cell;
}
template<typename T, unsigned int N>
__device__ T getVariable(const char(&variable_name)[N]) const;
template<typename T, flamegpu::size_type N, unsigned int M>
__device__ T getVariable(const char(&variable_name)[M], unsigned int index) const;
};
class iterator {
Message _message;
public:
__device__ iterator(const Filter& parent, const int relative_x)
: _message(parent, relative_x) {
// Increment to find first message
++_message;
}
__device__ iterator& operator++() { ++_message; return *this; }
__device__ iterator operator++(int) {
iterator temp = *this;
++* this;
return temp;
}
__device__ bool operator==(const iterator& rhs) const { return _message == rhs._message; }
__device__ bool operator!=(const iterator& rhs) const { return _message != rhs._message; }
__device__ Message& operator*() { return _message; }
__device__ Message* operator->() { return &_message; }
};
inline __device__ Filter(size_type _length, size_type x, size_type _radius);
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
inline __device__ Filter();
#endif
inline __device__ iterator begin(void) const {
// Bin before initial bin, as the constructor calls increment operator
return iterator(*this, min_cell - 1);
}
inline __device__ iterator end(void) const {
// Final bin, as the constructor calls increment operator
return iterator(*this, max_cell);
}
private:
size_type loc;
int min_cell;
int max_cell;
const size_type length;
};
__device__ In(const void *metadata)
: length(reinterpret_cast<const MetaData*>(metadata)->length)
{ }
inline __device__ WrapFilter wrap(const size_type x, const size_type radius = 1) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (radius == 0) {
DTHROW("Invalid radius %u for accessing array messagelist of length %u\n", radius, length);
return WrapFilter();
} else if ((radius * 2) + 1 > length) {
unsigned int min_r = length % 2 == 0 ? length - 2 : length - 1;
min_r /= 2;
DTHROW("%u is not a valid radius for accessing Array message lists, as the diameter of messages accessed exceeds the message list length (%u)."
" Maximum supported radius for this message list is %u.\n",
radius, length, min_r);
return WrapFilter();
} else if (x >= length) {
DTHROW("%u is not a valid position for iterating an Array message list of length %u, location must be within bounds.",
x, length);
return WrapFilter();
}
#endif
return WrapFilter(length, x, radius);
}
inline __device__ Filter operator() (const size_type x, const size_type radius = 1) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (radius == 0) {
DTHROW("Invalid radius %u for accessing array messagelist of length %u\n", radius, length);
return Filter();
} else if (x >= length) {
DTHROW("%u is not a valid position for iterating an Array message list of length %u, location must be within bounds.",
x, length);
return Filter();
}
#endif
return Filter(length, x, radius);
}
__device__ size_type size(void) const {
return length;
}
__device__ Message at(const size_type index) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (index >= length) {
DTHROW("Index is out of bounds for Array messagelist (%u >= %u).\n", index, length);
return Message(*this);
}
#endif
return Message(*this, index);
}
private:
const size_type length;
};
class MessageArray::Out {
public:
__device__ Out(const void *_metadata, unsigned int *scan_flag_messageOutput)
: scan_flag(scan_flag_messageOutput)
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
, metadata(reinterpret_cast<const MetaData*>(_metadata))
#else
, metadata(nullptr)
#endif
{ }
__device__ inline void setIndex(const size_type id) const;
template<typename T, unsigned int N>
__device__ void setVariable(const char(&variable_name)[N], T value) const;
template<typename T, unsigned int N, unsigned int M>
__device__ void setVariable(const char(&variable_name)[M], unsigned int index, T value) const;
protected:
unsigned int *scan_flag;
const MetaData * const metadata;
};
template<typename T, unsigned int N>
__device__ T MessageArray::In::Message::getVariable(const char(&variable_name)[N]) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
return detail::curve::DeviceCurve::getMessageVariable<T>(variable_name, index);
}
template<typename T, flamegpu::size_type N, unsigned int M> __device__
T MessageArray::In::Message::getVariable(const char(&variable_name)[M], const unsigned int array_index) const {
// simple indexing assumes index is the thread number (this may change later)
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
T value = detail::curve::DeviceCurve::getMessageArrayVariable<T, N>(variable_name, index, array_index);
return value;
}
template<typename T, unsigned int N>
__device__ T MessageArray::In::WrapFilter::Message::getVariable(const char(&variable_name)[N]) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index_1d >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
return detail::curve::DeviceCurve::getMessageVariable<T>(variable_name, index_1d);
}
template<typename T, flamegpu::size_type N, unsigned int M> __device__
T MessageArray::In::WrapFilter::Message::getVariable(const char(&variable_name)[M], const unsigned int array_index) const {
// simple indexing assumes index is the thread number (this may change later)
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index_1d >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
T value = detail::curve::DeviceCurve::getMessageArrayVariable<T, N>(variable_name, index_1d, array_index);
return value;
}
template<typename T, unsigned int N>
__device__ T MessageArray::In::Filter::Message::getVariable(const char(&variable_name)[N]) const {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index_1d >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
return detail::curve::DeviceCurve::getMessageVariable<T>(variable_name, index_1d);
}
template<typename T, flamegpu::size_type N, unsigned int M> __device__
T MessageArray::In::Filter::Message::getVariable(const char(&variable_name)[M], const unsigned int array_index) const {
// simple indexing assumes index is the thread number (this may change later)
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
// Ensure that the message is within bounds.
if (index_1d >= this->_parent.length) {
DTHROW("Invalid Array message, unable to get variable '%s'.\n", variable_name);
return static_cast<T>(0);
}
#endif
// get the value from curve using the message index.
T value = detail::curve::DeviceCurve::getMessageArrayVariable<T, N>(variable_name, index_1d, array_index);
return value;
}
template<typename T, unsigned int N>
__device__ void MessageArray::Out::setVariable(const char(&variable_name)[N], T value) const { // message name or variable name
if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in MessageArray::Out::setVariable().\n", variable_name);
#endif
return; // Fail silently
}
unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;
// set the variable using curve
detail::curve::DeviceCurve::setMessageVariable<T>(variable_name, value, index);
// setIndex() sets the optional message scan flag
}
template<typename T, unsigned int N, unsigned int M>
__device__ void MessageArray::Out::setVariable(const char(&variable_name)[M], const unsigned int array_index, T value) const {
if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in MessageArray::Out::setVariable().\n", variable_name);
#endif
return; // Fail silently
}
unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;
// set the variable using curve
detail::curve::DeviceCurve::setMessageArrayVariable<T, N>(variable_name, value, index, array_index);
// setIndex() sets the optional message scan flag
}
__device__ void MessageArray::Out::setIndex(const size_type id) const {
unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (id >= metadata->length) {
DTHROW("MessageArray index [%u] is out of bounds [%u]\n", id, metadata->length);
}
#endif
// set the variable using curve
detail::curve::DeviceCurve::setMessageVariable<size_type>("___INDEX", id, index);
// Set scan flag incase the message is optional
this->scan_flag[index] = 1;
}
__device__ MessageArray::In::WrapFilter::WrapFilter(const size_type _length, const size_type x, const size_type _radius)
: radius(_radius)
, length(_length) {
loc = x;
}
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
__device__ inline MessageArray::In::WrapFilter::WrapFilter()
: radius(0)
, length(0) {
loc = 0;
}
#endif
__device__ MessageArray::In::WrapFilter::Message& MessageArray::In::WrapFilter::Message::operator++() {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (!_parent.length)
return *this;
#endif
relative_cell++;
// Skip origin cell
if (relative_cell == 0) {
relative_cell++;
}
// Wrap over boundaries
index_1d = (this->_parent.loc + relative_cell + this->_parent.length) % this->_parent.length;
return *this;
}
__device__ MessageArray::In::Filter::Filter(const size_type _length, const size_type x, const size_type _radius)
: length(_length) {
loc = x;
min_cell = static_cast<int>(x) - static_cast<int>(_radius) < 0 ? -static_cast<int>(x) : -static_cast<int>(_radius);
max_cell = x + _radius >= _length ? static_cast<int>(_length) - 1 - static_cast<int>(x) : static_cast<int>(_radius);
}
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
__device__ inline MessageArray::In::Filter::Filter()
: length(0) {
loc = 0;
min_cell = 1;
max_cell = 0;
}
#endif
__device__ MessageArray::In::Filter::Message& MessageArray::In::Filter::Message::operator++() {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
if (!_parent.length)
return *this;
#endif
relative_cell++;
// Skip origin cell
if (relative_cell == 0) {
relative_cell++;
}
// Solve to 1 dimensional bin index
index_1d = this->_parent.loc + relative_cell;
return *this;
}
} // namespace flamegpu
#endif // INCLUDE_FLAMEGPU_RUNTIME_MESSAGING_MESSAGEARRAY_MESSAGEARRAYDEVICE_CUH_