Project
Loading...
Searching...
No Matches
InteractionSampler.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
14#ifndef ALICEO2_INTERACTIONSAMPLER_H
15#define ALICEO2_INTERACTIONSAMPLER_H
16
17#include <Rtypes.h>
18#include <TMath.h>
19#include <TRandom.h>
20#include <vector>
25
26namespace o2
27{
28namespace steer
29{
31{
32 public:
33 static constexpr float Sec2NanoSec = 1.e9; // s->ns conversion
35 void generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest);
36
37 void init();
38
39 void setInteractionRate(float rateHz)
40 {
41 mIntRate = rateHz;
42 mMuBC = -1.; // invalidate
43 }
44 float getInteractionRate() const { return mIntRate; }
46 {
47 mFirstIR.InteractionRecord::operator=(ir);
48 if (mFirstIR.orbit == 0 && mFirstIR.bc < 4) {
49 mFirstIR.bc = 4;
50 }
51 }
52 const o2::InteractionRecord& getFirstIR() const { return mFirstIR; }
53
54 void setMuPerBC(float mu)
55 {
56 mMuBC = mu;
57 mIntRate = -1.; // invalidate
58 }
59 float getMuPerBC() const { return mMuBC; }
60 void setBCTimeRMS(float tNS = 0.2) { mBCTimeRMS = tNS; }
61 float getBCTimeRMS() const { return mBCTimeRMS; }
62 const BunchFilling& getBunchFilling() const { return mBCFilling; }
65 void setBunchFilling(const std::string& bcFillingFile);
66
67 void print() const;
68
69 protected:
71 void nextCollidingBC(int n);
72
73 o2::math_utils::RandomRing<10000> mBCJumpGenerator; // generator of random jumps in BC
74 o2::math_utils::RandomRing<1000> mNCollBCGenerator; // generator of number of interactions in BC
75 o2::math_utils::RandomRing<1000> mCollTimeGenerator; // generator of number of interactions in BC
76
79 int mIntBCCache = 0;
80
81 float mIntRate = -1.;
82 float mBCTimeRMS = 0.2;
83 double mMuBC = -1.;
84
86 std::vector<float> mTimeInBC;
87 std::vector<uint16_t> mInteractingBCs; // vector of interacting BCs
88 int mCurrBCIdx = 0;
89
90 static constexpr float DefIntRate = 50e3;
91
93};
94
95//_________________________________________________
96inline void InteractionSampler::generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest)
97{
98 // fill vector with interaction records
99 dest.clear();
100 for (int i = dest.capacity(); i--;) {
101 dest.push_back(generateCollisionTime());
102 }
103}
104
105//_________________________________________________
107{
109 if ((mCurrBCIdx += n) >= (int)mInteractingBCs.size()) {
111 mCurrBCIdx %= mInteractingBCs.size();
112 }
114}
115
116} // namespace steer
117} // namespace o2
118
119#endif
uint64_t bc
Definition RawEventData.h:5
int32_t i
Header to collect LHC related constants.
o2::InteractionTimeRecord mFirstIR
static constexpr float DefIntRate
default interaction rate
static constexpr float Sec2NanoSec
std::vector< uint16_t > mInteractingBCs
o2::BunchFilling mBCFilling
patter of active BCs
o2::math_utils::RandomRing< 10000 > mBCJumpGenerator
float mBCTimeRMS
BC time spread in NANOSECONDS.
double mMuBC
interaction probability per BC
const BunchFilling & getBunchFilling() const
o2::math_utils::RandomRing< 1000 > mNCollBCGenerator
o2::math_utils::RandomRing< 1000 > mCollTimeGenerator
void generateCollisionTimes(std::vector< o2::InteractionTimeRecord > &dest)
std::vector< float > mTimeInBC
interaction times within single BC
int mCurrBCIdx
counter for current interacting bunch
void setFirstIR(const o2::InteractionRecord &ir)
const o2::InteractionTimeRecord & generateCollisionTime()
o2::InteractionTimeRecord mIR
const o2::InteractionRecord & getFirstIR() const
void setBunchFilling(const BunchFilling &bc)
ClassDefNV(InteractionSampler, 1)
float mIntRate
total interaction rate in Hz
int mIntBCCache
N interactions left for current BC.
void setInteractionRate(float rateHz)
GLdouble n
Definition glcorearb.h:1982
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
uint32_t orbit
LHC orbit.
uint16_t bc
bunch crossing ID of interaction
o2::InteractionRecord ir(0, 0)