EICd
EIC data model
EICd Documentation

EICd - EIC data model

Overview

A podio based data model.

Full Description File

The entire data model is defined with a single YAML file. Here is the current definition:

---
options :
  # should getters / setters be prefixed with get / set?
  getSyntax: False
  # should POD members be exposed with getters/setters in classes that have them as members?
  exposePODMembers: False
  includeSubfolder: True

## right now we rigurously enforce:
##  - No breaking of PODness:
##        - No use of relations and vectors
##     - Use special Relation structures where needed, indexing by ID (index)
##        - IDs are stored as eic::Index, which is a thin layer of an signed integer
##          where -1 relates to "no index".
##        - For 1-to-many relations: Use many-to-1 IDs instead --> use forward links
##          This puts the burden on the reconstruction algorithm and keeps the data 2D!
##  - Use float most of the time except for 4-vectors where ppm precision is important.
##  - Data alignment: 
##        - data should be aligned with a 64-bit structure where possible.
##        - when using 32 bit values, use them in pairs (or after all 64-bit variables are defined). 
##        - same goes for 16-bit values (keep them aligned with the largest following component)
##  - Explicitly specify the integer length (use the typedefs from <cstdint>, 
##    such as int32_t etc)

components:

  ## Unique field identifier. Has 2 components: source and ID where
  ## source identifies the originating collection (or algorithm) and ID 
  ## the ID of the entry within this collection. 
  ## Defaults to -1 for an unset index.
  eic::Index:
    Members:
      - int32_t    value
      - int32_t    source
    ExtraCode:
      declaration: "
        Index() : source{0}, value{-1} {}\n
        Index(int32_t id, int32_t src) : value{id}, source{src} {}\n
        Index(Index rhs, int32_t new_src) : Index(rhs.value, new_src) {} \n
        bool empty() const {return value < 0;}\n
        bool valid() const {return value >= 0 && source != 0;}\n
        bool equals(const eic::Index& rhs) const {return rhs.source == source && rhs.value == value;}\n
        bool operator==(const eic::Index& rhs) const {return equals(rhs);}\n
        bool operator!=(const eic::Index& rhs) const {return !equals(rhs);}\n
        explicit operator bool() const {return valid();}
      "

  ## simple weight that defaults to 1 if not set
  eic::Weight:
    Members:
      - float value
    ExtraCode:
      declaration: "
        Weight() : value{1.} {}\n
        Weight(double w) : value {static_cast<float>(w)} {}\n
        Weight& operator=(double w) {value = static_cast<float>(w); return *this;}\n
        operator float() const {return value;}
      "

  ## first-second pair of float s
  eic::FloatPair:
    Members:
      - float first
      - float second
    ExtraCode:
      includes: "#include <tuple>"
      declaration: "
        FloatPair() : first{0}, second{0} {}\n
        FloatPair(double  a, double  b) : first{static_cast<float>(a)}, second{static_cast<float>(b)} {}\n
        FloatPair(const std::pair<float, float>& p) : first{p.first}, second{p.second} {}\n
        FloatPair& operator=(const std::pair<float, float>& p) {first = p.first; second = p.second; return *this;}\n
        float& operator[](unsigned i) {return *(&first + i);}\n
        const float& operator[](unsigned i) const {return *(&first + i);}\n
        operator std::pair<float, float>() const {return {first, second};}\n
      "

  eic::VectorXY:
    Members:
      - float x          // [mm] or [GeV]
      - float y          //
    ExtraCode:
      includes: "#include <cmath>\n"
      declaration: "
        VectorXY() : x{0}, y{0} {}\n
        VectorXY(double xx, double yy) : x{static_cast<float>(xx)}, y{static_cast<float>(yy)} {}\n
        float& operator[](unsigned i) {return *(&x + i);}\n
        const float& operator[](unsigned i) const {return *(&x + i);}\n
        float mag() const {return std::hypot(x, y);}\n
        float r() const {return mag();}\n
        float phi() const {return atan2(y, x);}\n
        operator std::pair<float, float>() const {return {x, y};}\n
        float dot(const VectorXY& rhs) const {return x*rhs.x + y*rhs.y;}\n
        VectorXY add(const VectorXY& rhs) const {return {x+rhs.x, y+rhs.y};}\n
        VectorXY subtract(const VectorXY& rhs) const {return {x-rhs.x, y-rhs.y};}\n
        VectorXY scale(double f) const {return {f*x, f*y};}\n
      "

  eic::Direction:
    Members:
      - float theta      // [rad, 0->pi]
      - float phi        // [rad, -pi->pi]
    ExtraCode:
      includes: "#include <cmath>\n#include <tuple>"
      declaration: "
        Direction() : theta{0}, phi{0} {}\n
        Direction(double th, double ph) : theta{static_cast<float>(th)}, phi{static_cast<float>(ph)} {}\n
        Direction(double x, double y, double z)\n
          : theta{static_cast<float>(acos(z/std::hypot(x,y,z)))}\n
          , phi{static_cast<float>(atan2(y,x))} {}\n
        template <class VectorType> Direction(const VectorType& v) : Direction(v.theta(), v.phi()) {}\n
        operator std::pair<float, float>() const {return {theta, phi};}\n
        float eta() const {return -log(tan(0.5*theta));}\n
        Direction add(const Direction& rhs) const {return {theta+rhs.theta, phi+rhs.phi};}\n
        Direction subtract(const Direction& rhs) const {return {theta-rhs.theta, phi-rhs.phi};}\n
      "

  eic::VectorXYZ:
    Members:
      - float x          // [mm] or [GeV]
      - float y          //
      - float z          //
    ExtraCode:
      includes: "#include <cmath>\n#include<tuple>"
      declaration: "
        VectorXYZ() : x{0}, y{0}, z{0} {}\n
        VectorXYZ(double xx, double yy, double zz) : x{static_cast<float>(xx)}, y{static_cast<float>(yy)}, z{static_cast<float>(zz)} {}\n
        template<class VectorPolarType> VectorXYZ(const VectorPolarType& v) : VectorXYZ(v.x(), v.y(), v.z()) {}\n
        float& operator[](unsigned i) {return *(&x + i);}\n
        const float& operator[](unsigned i) const {return *(&x + i);}\n
        float mag() const {return std::hypot(x, y, z);}\n
        float r() const {return mag();}\n
        float theta() const {return acos(z/mag());}\n
        float phi() const {return atan2(y,x);}\n
        float eta() const {return -log(tan(0.5*theta()));}\n
        operator std::tuple<float, float, float>() {return {x, y, z};}\n
        float dot(const VectorXYZ& rhs) const {return x*rhs.x + y*rhs.y + z*rhs.z;}\n
        VectorXYZ add(const VectorXYZ& rhs) const {return {x+rhs.x, y+rhs.y, z+rhs.z};}\n
        VectorXYZ subtract(const VectorXYZ& rhs) const {return {x-rhs.x, y-rhs.y, z-rhs.z};}\n
        VectorXYZ scale(double f) const {return {f*x, f*y, f*z};}\n
      "
  eic::VectorPolar:
    Members:
      - float r          // [mm] or [GeV]
      - float theta      // [rad, 0->pi]
      - float phi        // [rad, -pi->pi]
    ExtraCode:
      includes: "#include <cmath>\n#include<tuple>"
      declaration: "
        VectorPolar() : r{0}, theta{0}, phi{0} {}\n
        VectorPolar(double rr, double th, double ph) : r{static_cast<float>(rr)}, theta{static_cast<float>(th)}, phi{static_cast<float>(ph)} {}\n
        template<class VectorXYZType> VectorPolar(const VectorXYZType& v) : VectorPolar(v.r(), v.theta(), v.phi()) {}\n
        float mag() const {return r;}\n
        float x() const {return r * cos(phi) * sin(theta);}\n
        float y() const {return r * sin(phi) * sin(theta);}\n
        float z() const {return r * cos(theta);}\n
        float eta() const {return -log(tan(0.5*theta));}\n
        operator std::tuple<float, float, float>() {return {r, theta, phi};}\n
      "

  eic::VectorXYZT:
    Members:
      - double x        // [mm] or [GeV]
      - double y        //
      - double z        //
      - double t        // [ns] or [GeV]
    ExtraCode:
      includes: "#include <cmath>\n#include <tuple>"
      declaration: "
        VectorXYZT() : x{0}, y{0}, z{0}, t{0} {}\n
        VectorXYZT(double xx, double yy, double zz, double tt) : x{xx}, y{yy}, z{zz}, t{tt} {}\n
        double& operator[](unsigned i) {return *(&x + i);}\n
        const double& operator[](unsigned i) const {return *(&x + i);}\n
        double mag() const {return std::hypot(x, y, z);}\n
        double r() const {return mag();}\n
        double theta() const {return acos(z/mag());}\n
        double phi() const {return atan2(y,x);}\n
        double eta() const {return -log(tan(0.5*theta()));}
        double energy() const {return t;}\n
        double mass() const {return sqrt(t*t - x*x - y*y - z*z);}\n
        operator std::tuple<double, double, double, double>() {return {x, y, z, t};}\n
        double dot(const VectorXYZT& rhs) const {return t*rhs.t - x*rhs.x - y*rhs.y - z*rhs.z;}\n
        VectorXYZT add(const VectorXYZT& rhs) const {return {x+rhs.x, y+rhs.y, z+rhs.z, t+rhs.t};}\n
        VectorXYZT subtract(const VectorXYZT& rhs) const {return {x-rhs.x, y-rhs.y, z-rhs.z, t-rhs.t};}\n
        VectorXYZT scale(double f) const {return {f*x, f*y, f*z, f*t};}\n
      "

  eic::CovDiagXYZ:
    Members:
      - float xx
      - float yy
      - float zz
    ExtraCode:
      declaration: "
        CovDiagXYZ() : xx{0}, yy{0}, zz{0} {}\n
        CovDiagXYZ(double x, double y, double z) : xx{static_cast<float>(x)}, yy{static_cast<float>(y)}, zz{static_cast<float>(z)} {}\n
        float operator()(unsigned i, unsigned j) const {return (i == j) ? *(&xx + i) : 0.;}\n
        "

  eic::CovDiagXYZT:
    Members:
      - double xx
      - double yy
      - double zz
      - double tt
    ExtraCode:
      declaration: "
        CovDiagXYZT() : xx{0}, yy{0}, zz{0}, tt{0} {}\n
        CovDiagXYZT(double x, double y, double z, double t) : xx{x}, yy{y}, zz{z}, tt{t} {}\n
        double operator()(unsigned i, unsigned j) const {return (i == j) ? *(&xx + i) : 0.;}\n
        "

  eic::CovXYZ:
    Members:
      - float xx
      - float yy
      - float zz
      - float xy
      - float xz
      - float yz
    ExtraCode:
      declaration: "
        CovXYZ() : xx{0}, yy{0}, zz{0}, xy{0}, xz{0}, yz{0} {}\n
        CovXYZ(double vx, double vy, double vz, double vxy, double vxz, double vyz)\n
          : xx{static_cast<float>(vx)}, yy{static_cast<float>(vy)}, zz{static_cast<float>(vz)},\n
            xy{static_cast<float>(vxy)}, xz{static_cast<float>(vxz)}, yz{static_cast<float>(vyz)} {}\n
        float operator()(unsigned i, unsigned j) const {\n
          // diagonal\n
          if (i == j) {\n
            return *(&xx + i);\n
          }\n
          // off-diagonal\n
          // we have as options (0, 1), (0, 2) and (1, 2) (and mirrored)\n
          // note that, starting from xy, we find the correct element at (i+j-1)\n
          return *(&xy + i + j - 1);\n
        }\n
      "

datatypes:

  ## ==========================================================================
  ## Event info
  ## ==========================================================================

  eic::EventInfo:
    Description: "Event Info"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - uint64_t          run              // Run number.
      - uint64_t          number           // Event number.
      - int32_t           type             // event type identifier (TBD).
      - int32_t           proc             // Process identifier (TBD).
      - int32_t           source           // Source/identifier (TBD)
      - eic::Weight       weight           // Optional event weight (useful for MC)

  
  ## ==========================================================================
  ## Particle info
  ## ==========================================================================

  eic::BasicParticle:
    Description: "Basic particle used internally to communicate basic particle properties."
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // Unique particle index
      - eic::VectorXYZ    p                 // momentum [GeV]
      - eic::VectorXYZ    v                 // vertex [mm]
      - float             time              // Time in [ns]
      - int32_t           pid               // particle PDG code
      - int16_t           status            // Status code
      - int16_t           charge            // Particle charge (or sign)
      - eic::Weight       weight            // Particle weight, e.g. from PID algorithm [0-1]

  eic::ReconstructedParticle:
    Description: "EIC Reconstructed Particle"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // Unique particle index
      - eic::VectorXYZ    p                 // momentum vector [GeV]
      - eic::VectorXYZ    v                 // vertex [mm]
      - float             time              // Time in [ns]
      - int32_t           pid               // PID of reconstructed particle.
      - int16_t           status            // Status code
      - int16_t           charge            // Particle charge (or sign)
      - eic::Weight       weight            // Particle weight, e.g. from PID algorithm [0-1]
      - eic::Direction    direction         // Direction (theta/phi of this particle [mrad])
      - float             momentum          // particle 3-momentum magnitude [GeV]
      - float             energy            // Particle energy, consistent with PID assigment [GeV]
      - float             mass              // The mass of the particle in [GeV]

  eic::ReconstructedParticleRelations:
    Description: "Relational info associated with our reconstructed particle"
    Author: "S. Joosten"
    Members:
      - eic::Index        recID             // ReconstructedParticle index
      - eic::Index        vertexID          // Start vertex for this particle
      - eic::Index        trackID           // Index of the associated track, if any
      - eic::Index        ecalID            // Index of associated pos/barrel/neg ECAL cluster, if any
      - eic::Index        hcalID            // Index of associated pos/barrel/neg HCAL cluster, if any
      - eic::Index        cherID            // Index of associated pos/barrel/neg Cherenkov info, if any
      - eic::Index        tofID             // Index of the associated TOF info, if any
      - eic::Index        mcID              // Index of the associated MC particle, if any

  ## ==========================================================================
  ## Calorimetry
  ## ==========================================================================
  eic::RawCalorimeterHit:
    Description: "Raw (digitized) calorimeter hit"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique ID for this hit
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - int64_t           amplitude         // The amplitude of the hit in ADC counts.
      - int64_t           time              // Timing in TDC

  eic::CalorimeterHit:
    Description: "Calorimeter hit"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique ID for this hit
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - int32_t           layer             // layer for this hit
      - int32_t           sector            // sector for this hit
      - float             energy            // The energy for this hit in [GeV].
      - float             energyError       // Error on energy [GeV].
      - float             time              // The time of the hit in [ns].
      - eic::VectorXYZ    position          // The global position of the hit in world coordinates [mm].
      - eic::VectorXYZ    local             // The local position of the hit in detector coordinates [mm].
      - eic::VectorXYZ    dimension         // The dimension information of the cell [mm].


  ## ==========================================================================
  ## Clustering
  ## ==========================================================================
  
  eic::ProtoCluster:
    Description: "Relational info linking hits to their associated cluster"
    Author: "S. Joosten"
    Members:
      - eic::Index        hitID             // Hit ID
      - eic::Index        clusterID         // ID of the cluster associated with this hit (-1 if none)
      - eic::Weight       weight            // How much of this hit belongs to the cluster? [0->1]

  eic::Cluster:
    Description: "EIC cluster"
    Author: "W. Armstrong, S. Joosten, C.Peng"
    Members:
      - eic::Index        ID                // unique ID for this cluster
      - float             energy            // Reconstructed energy of the cluster [GeV].
      - float             energyError       // Error on the cluster energy [GeV]
      - float             time              // [ns]
      - uint32_t          nhits             // Number of hits in the cluster.
      - eic::VectorXYZ    position          // Global position of the cluster [mm].
      - eic::CovXYZ       positionError     // Covariance matrix of the position (6 Parameters).
      - float             radius            // shower radius [mm]
      - float             skewness          // shower skewness [unitless]

  eic::Cluster2DInfo:
    Description: "Additional info for 3D clusters" 
    Author: "S. Joosten"
    Members:
      - eic::Index        clusterID         // Primary cluster ID
      - eic::VectorPolar  polar             // Cluster position polar information
      - float             eta               // Cluster pseudorapidity

  eic::Cluster3DInfo:
    Description: "Additional info for 3D clusters" 
    Author: "S. Joosten"
    Members:
      - eic::Index        clusterID         // Primary cluster ID
      - eic::VectorPolar  polar             // Cluster position polar information
      - float             eta               // Cluster pseudorapidity
      - eic::Direction    direction         // Intrinsic direction of the cluster at the central position [rad, 0->pi and -pi->pi]

  eic::ClusterLayer:
    Description: "2D Cluster in a single layer for a multi-layer detector"
    Author: "S. Joosten, C. Peng"
    Members:
      - eic::Index        ID                // unique layer ID
      - eic::Index        clusterID         // associated full 3D cluster, -1 if none
      - int32_t           layer             // layer number for this cluster layer
      - uint32_t          nhits             // Number of hits
      - float             energy            // Energy in this cluster layer [GeV]
      - float             energyError       // Error on energy [GeV]
      - float             radius            // Shower radius [mm]
      - float             skewness          // Skewness of hits distribution
      - eic::VectorXYZ    position          // Global center position. [mm]

  eic::MergedClusterRelations:
    Description: "Relational info between a merged cluster and its parents"
    Author: "S. Joosten"
    Members:
      - eic::Index        clusterID         // associated cluster ID
      - uint32_t          size              // number of valid parents
      - std::array<eic::Index, 4> parent    // (up to 4) parents for this cluster

  ## ==========================================================================
  ## RICH/Cherenkov data structures
  ## ==========================================================================

  eic::RawPMTHit:
    Description: "EIC Raw PMT hit"
    Author: "S. Joosten, C. Peng"
    Members:
      - eic::Index        ID                // unique hit ID
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - uint32_t          amplitude         // PMT signal amplitude [ADC]
      - uint32_t          time              // PMT signal time [TDC]

  eic::PMTHit:
    Description: "EIC PMT hit"
    Author: "S. Joosten, C. Peng"
    Members:
      - eic::Index        ID                // Unique hit ID
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - float             npe               // estimated number of photo-electrons [#]
      - float             time              // Time [ns]
      - float             timeError         // Error on the time [ns]
      - eic::VectorXYZ    position          // PMT hit position [mm]
      - eic::VectorXYZ    local             // The local position of the hit in detector coordinates [mm]
      - eic::VectorXYZ    dimension         // The dimension information of the pixel [mm].

  eic::RingImage:
    Description: "EIC Ring Image Cluster"
    Author: "S. Joosten, C. Peng"
    Members:
      - eic::Index        ID                // Unique cluster ID
      - float             npe               // number of photo-electrons [#]
      - eic::VectorXYZ    position          // Global position of the cluster [mm]
      - eic::VectorXYZ    positionError     // Error on the position
      - float             theta             // opening angle of the ring [rad, 0->pi]
      - float             thetaError        // error on the opening angle
      - float             radius            // radius of the best fit ring [mm]
      - float             radiusError       // estimated error from the fit [mm]

  ## ==========================================================================
  ## Tracking
  ## ==========================================================================
  
  eic::RawTrackerHit:
    Description: "Raw (digitized) tracker hit"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique ID for this hit
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - int32_t           time              // tdc value.
      - int32_t           charge            // adc value

  eic::TrackerHit:
    Description: "Tracker hit (reconstructed from Raw)"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique ID for this hit
      - int64_t           cellID            // The detector specific (geometrical) cell id.
      - eic::VectorXYZT   position          // Hit (cell) position and time [mm, ns]
      - eic::CovDiagXYZT  covMatrix         // Covariance Matrix
      - float             edep              // Energy deposit in this hit [GeV]
      - float             edepError         // Error on the energy deposit [GeV]
    ConstExtraCode:
      declaration: "
        double time() const {return position().t;}\n
      "

    ## Here's a stub for a prototrack setup. If this is all we use we should
    ## probably just use protocluster instead, but I assume there will be other
    ## members we'll need to communicate
    ## eic::ProtoTrack:
    ##   Description: "Proto track info"
    ##   Author: "S. Joosten"
    ##   Members:
    ##     eic::Index       hitID             // Unique hit identifier
    ##     eic::Index       trackID           // link to the associated track
    ##     eic::Weight      weight            // prototrack weight, in case we share pixels [0-1]
  
  eic::TrackParameters:
    Description: "ACTS Bound Track parameters"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique ID for this track 
      - eic::FloatPair    loc               // tracking location
      - eic::FloatPair    locError          // error on the location
      - eic::Direction    direction         // track direction (theta, phi) [rad, 0-pi and -pi->pi]
      - eic::Direction    directionError    // error on the direction [rad]
      - float             qOverP            // [e/GeV]
      - float             qOverPError       // error on qOverP
      - float             time              // track time [ns]    
      - float             timeError         // error on the time

  ## ==========================================================================
  ## Vertexing
  ## ==========================================================================

  eic::Vertex:
    Description: "EIC vertex"
    Author: "W. Armstrong, S. Joosten"
    Members:
      - eic::Index        ID                // unique vertex ID
      - eic::VectorXYZ    position          // postion of vertex [mm]
      - float             time              // time of vertex [ns]
      - float             chi2              // Chi squared of the vertex fit.
      - float             probability       // Probability of the vertex fit
      - bool              primary           // Whether it is the primary vertex of the event

Installing

mkdir build && cd build
cmake ../. -DCMAKE_INSTALL_PREFIX=$HOME/stow/eicd -DBUILD_DATA_MODEL=ON
make -j4 install