Program Listing for File DeviceAPI.cuh

Return to documentation for file (include/flamegpu/runtime/DeviceAPI.cuh)

#ifndef INCLUDE_FLAMEGPU_RUNTIME_DEVICEAPI_CUH_
#define INCLUDE_FLAMEGPU_RUNTIME_DEVICEAPI_CUH_


#include <cassert>
#include <cstdint>
#include <limits>

#ifndef __CUDACC_RTC__
#include "flamegpu/runtime/detail/curve/DeviceCurve.cuh"
#include "flamegpu/runtime/messaging_device.h"
#else
#include "dynamic/curve_rtc_dynamic.h"
#endif  // !_RTC
#include "flamegpu/runtime/random/AgentRandom.cuh"
#include "flamegpu/runtime/environment/DeviceEnvironment.cuh"
#include "flamegpu/runtime/AgentFunction.cuh"
#include "flamegpu/runtime/AgentFunctionCondition.cuh"
#include "flamegpu/defines.h"

#ifdef FLAMEGPU_USE_GLM
#ifdef __CUDACC__
#ifdef __NVCC_DIAG_PRAGMA_SUPPORT__
#pragma nv_diag_suppress = esa_on_defaulted_function_ignored
#else
#pragma diag_suppress = esa_on_defaulted_function_ignored
#endif  // __NVCC_DIAG_PRAGMA_SUPPORT__
#endif  // __CUDACC__
#include <glm/glm.hpp>
#endif  // FLAMEGPU_USE_GLM

namespace flamegpu {

class ReadOnlyDeviceAPI {
    // Friends have access to TID() & TS_ID()
    template<typename AgentFunctionCondition>
    friend __global__ void agent_function_condition_wrapper(
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
        exception::DeviceExceptionBuffer *,
#endif
#ifndef __CUDACC_RTC__
        const detail::curve::CurveTable *,
#endif
        const unsigned int,
        detail::curandState *,
        unsigned int *);

 public:
    __device__ ReadOnlyDeviceAPI(detail::curandState *&d_rng)
        : random(AgentRandom(&d_rng[getIndex()]))
        , environment(DeviceEnvironment()) { }
    template<typename T, unsigned int N> __device__
    T getVariable(const char(&variable_name)[N]) const;
    template<typename T, unsigned int N, unsigned int M> __device__
    T getVariable(const char(&variable_name)[M], unsigned int index) const;
    __device__ id_t getID() {
        return getVariable<id_t>("_id");
    }

    __forceinline__ __device__ unsigned int getStepCounter() const {
        return environment.getProperty<unsigned int>("_stepCount");
    }

    const AgentRandom random;
    const ReadOnlyDeviceEnvironment environment;

    __forceinline__ __device__ static unsigned int getIndex() {
        /*
        // 3D version
        auto blockId = blockIdx.x + blockIdx.y * gridDim.x
        + gridDim.x * gridDim.y * blockIdx.z;
        auto threadId = blockId * (blockDim.x * blockDim.y * blockDim.z)
        + (threadIdx.z * (blockDim.x * blockDim.y))
        + (threadIdx.y * blockDim.x)
        + threadIdx.x;
        return threadId;*/
#ifdef FLAMEGPU_SEATBELTS
        assert(blockDim.y == 1);
        assert(blockDim.z == 1);
        assert(gridDim.y == 1);
        assert(gridDim.z == 1);
#endif
        return blockIdx.x * blockDim.x + threadIdx.x;
    }
    __forceinline__ __device__ bool isAgent(const char* agent_name) {
        return detail::curve::DeviceCurve::isAgent(agent_name);
    }
    __forceinline__ __device__ bool isState(const char* agent_state) {
        return detail::curve::DeviceCurve::isState(agent_state);
    }
};

template<typename MessageIn, typename MessageOut>
class DeviceAPI {
    // Friends have access to TID() & TS_ID()
    template<typename AgentFunction, typename _MessageIn, typename _MessageOut>
    friend __global__ void agent_function_wrapper(
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
        exception::DeviceExceptionBuffer *,
#endif
#ifndef __CUDACC_RTC__
        const detail::curve::CurveTable *,
#endif
        id_t*,
        const unsigned int,
        const void *,
        const void *,
        detail::curandState *,
        unsigned int *,
        unsigned int *,
        unsigned int *);

 public:
    class AgentOut {
     public:
        __device__ AgentOut(id_t *&d_agent_output_nextID, unsigned int *&scan_flag_agentOutput)
            : scan_flag(scan_flag_agentOutput)
            , nextID(d_agent_output_nextID) { }
        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;
        __device__ id_t getID() const;

     private:
        __device__ void genID() const;
        unsigned int* const scan_flag;
        mutable id_t id = ID_NOT_SET;
        id_t *nextID;
    };
    __device__ DeviceAPI(
        id_t *&d_agent_output_nextID,
        detail::curandState *&d_rng,
        unsigned int *&scanFlag_agentOutput,
        typename MessageIn::In &&message_in,
        typename MessageOut::Out &&message_out)
        : message_in(message_in)
        , message_out(message_out)
        , agent_out(AgentOut(d_agent_output_nextID, scanFlag_agentOutput))
        , random(AgentRandom(&d_rng[getIndex()]))
        , environment(DeviceEnvironment())
    { }
    template<typename T, unsigned int N> __device__
    T getVariable(const char(&variable_name)[N]) const;
    template<typename T, unsigned int N, unsigned int M> __device__
    T getVariable(const char(&variable_name)[M], unsigned int index) const;
    template<typename T, unsigned int N>
    __device__ void setVariable(const char(&variable_name)[N], T value);
    template<typename T, unsigned int N, unsigned int M>
    __device__ void setVariable(const char(&variable_name)[M], unsigned int index, T value);
    __device__ id_t getID() {
        return getVariable<id_t>("_id");
    }

    __forceinline__ __device__ unsigned int getStepCounter() const {
        return environment.getProperty<unsigned int>("_stepCount");
    }

    __forceinline__ __device__ static unsigned int getIndex() {
        /*
        // 3D version
        auto blockId = blockIdx.x + blockIdx.y * gridDim.x
        + gridDim.x * gridDim.y * blockIdx.z;
        auto threadId = blockId * (blockDim.x * blockDim.y * blockDim.z)
        + (threadIdx.z * (blockDim.x * blockDim.y))
        + (threadIdx.y * blockDim.x)
        + threadIdx.x;
        return threadId;*/
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
        assert(blockDim.y == 1);
        assert(blockDim.z == 1);
        assert(gridDim.y == 1);
        assert(gridDim.z == 1);
#endif
        return blockIdx.x * blockDim.x + threadIdx.x;
    }
    __forceinline__ __device__ bool isAgent(const char* agent_name) {
        return detail::curve::DeviceCurve::isAgent(agent_name);
    }
    __forceinline__ __device__ bool isState(const char* agent_state) {
        return detail::curve::DeviceCurve::isState(agent_state);
    }

    const typename MessageIn::In message_in;
    const typename MessageOut::Out message_out;
    const AgentOut agent_out;
    const AgentRandom random;
    const DeviceEnvironment environment;
};


/******************************************************************************************************* Implementation ********************************************************/

template<typename T, unsigned int N>
__device__ T ReadOnlyDeviceAPI::getVariable(const char(&variable_name)[N]) const {
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

    // get the value from curve
    T value = detail::curve::DeviceCurve::getAgentVariable<T>(variable_name, index);

    // return the variable from curve
    return value;
}
template<typename T, unsigned int N, unsigned int M>
__device__ T ReadOnlyDeviceAPI::getVariable(const char(&variable_name)[M], const unsigned int array_index) const {
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

    // get the value from curve
    T value = detail::curve::DeviceCurve::getAgentArrayVariable<T, N>(variable_name, index, array_index);

    // return the variable from curve
    return value;
}

template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N>
__device__ T DeviceAPI<MessageIn, MessageOut>::getVariable(const char(&variable_name)[N]) const {
    using detail::sm;
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

    // get the value from curve
    T value = detail::curve::DeviceCurve::getAgentVariable<T>(variable_name, index);

    // return the variable from curve
    return value;
}

template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N, unsigned int M>
__device__ T DeviceAPI<MessageIn, MessageOut>::getVariable(const char(&variable_name)[M], const unsigned int array_index) const {
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

    // get the value from curve
    T value = detail::curve::DeviceCurve::getAgentArrayVariable<T, N>(variable_name, index, array_index);

    // return the variable from curve
    return value;
}

template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N>
__device__ void DeviceAPI<MessageIn, MessageOut>::setVariable(const char(&variable_name)[N], T value) {
    if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
        DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in DeviceAPI::setVariable().\n", variable_name);
#endif
        return;  // Fail silently
    }
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;
    // set the variable using curve
    detail::curve::DeviceCurve::setAgentVariable<T>(variable_name, value, index);
}
template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N, unsigned int M>
__device__ void DeviceAPI<MessageIn, MessageOut>::setVariable(const char(&variable_name)[M], const unsigned int array_index, const T value) {
    if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
        DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in DeviceAPI::setVariable().\n", variable_name);
#endif
        return;  // Fail silently
    }
    // simple indexing assumes index is the thread number (this may change later)
    const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

    // set the variable using curve
    detail::curve::DeviceCurve::setAgentArrayVariable<T, N>(variable_name, value, index, array_index);
}

template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N>
__device__ void DeviceAPI<MessageIn, MessageOut>::AgentOut::setVariable(const char(&variable_name)[N], T value) const {
    if (nextID) {
        if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
            DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in AgentOut::setVariable().\n", variable_name);
#endif
            return;  // Fail silently
        }
        // simple indexing assumes index is the thread number (this may change later)
        const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

        // set the variable using curve
        detail::curve::DeviceCurve::setNewAgentVariable<T>(variable_name, value, index);

        // Mark scan flag
        genID();
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
    } else {
        DTHROW("Agent output must be enabled per agent function when defining the model.\n");
#endif
    }
}
template<typename MessageIn, typename MessageOut>
template<typename T, unsigned int N, unsigned int M>
__device__ void DeviceAPI<MessageIn, MessageOut>::AgentOut::setVariable(const char(&variable_name)[M], const unsigned int array_index, T value) const {
    if (nextID) {
        if (variable_name[0] == '_') {
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
            DTHROW("Variable names starting with '_' are reserved for internal use, with '%s', in AgentOut::setVariable().\n", variable_name);
#endif
            return;  // Fail silently
        }
        // simple indexing assumes index is the thread number (this may change later)
        const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;

        // set the variable using curve
        detail::curve::DeviceCurve::setNewAgentArrayVariable<T, N>(variable_name, value, index, array_index);

        // Mark scan flag
        genID();
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
    } else {
        DTHROW("Agent output must be enabled per agent function when defining the model.\n");
#endif
    }
}

template<typename MessageIn, typename MessageOut>
__device__ id_t DeviceAPI<MessageIn, MessageOut>::AgentOut::getID() const {
    if (nextID) {
        genID();
        return this->id;
    }
#if !defined(FLAMEGPU_SEATBELTS) || FLAMEGPU_SEATBELTS
    DTHROW("Agent output must be enabled per agent function when defining the model.\n");
#endif
    return ID_NOT_SET;
}
#ifdef __CUDACC__
template<typename MessageIn, typename MessageOut>
__device__ void DeviceAPI<MessageIn, MessageOut>::AgentOut::genID() const {
    // Only called internally, so no need to check nextID != nullptr
    // Only assign id and scan flag once
    if (this->id == ID_NOT_SET) {
        this->id = atomicInc(this->nextID, std::numeric_limits<id_t>().max());
        const unsigned int index = (blockDim.x * blockIdx.x) + threadIdx.x;
        detail::curve::DeviceCurve::setNewAgentVariable<id_t>("_id", this->id, index);  // Can't use ID_VARIABLE_NAME inline, as it isn't of char[N] type
        this->scan_flag[index] = 1;
    }
}
#endif

}  // namespace flamegpu

#endif  // INCLUDE_FLAMEGPU_RUNTIME_DEVICEAPI_CUH_