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#include <TH1F.h>
26
27namespace o2
28{
29namespace steer
30{
32{
33 public:
34 static constexpr float Sec2NanoSec = 1.e9; // s->ns conversion
36 void generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest);
37
38 void init();
39
40 void setInteractionRate(float rateHz)
41 {
42 mIntRate = rateHz;
43 mMuBC = -1.; // invalidate
44 }
45 float getInteractionRate() const { return mIntRate; }
47 {
48 mFirstIR.InteractionRecord::operator=(ir);
49 if (mFirstIR.orbit == 0 && mFirstIR.bc < 4) {
50 mFirstIR.bc = 4;
51 }
52 }
53 const o2::InteractionRecord& getFirstIR() const { return mFirstIR; }
54
55 void setMuPerBC(float mu)
56 {
57 mMuBC = mu;
58 mIntRate = -1.; // invalidate
59 }
60 float getMuPerBC() const { return mMuBC; }
61 void setBCTimeRMS(float tNS = 0.2) { mBCTimeRMS = tNS; }
62 float getBCTimeRMS() const { return mBCTimeRMS; }
63 const BunchFilling& getBunchFilling() const { return mBCFilling; }
66 void setBunchFilling(const std::string& bcFillingFile);
67
68 void print() const;
69
70 protected:
71 virtual int simulateInteractingBC();
72 void nextCollidingBC(int n);
73
74 o2::math_utils::RandomRing<10000> mBCJumpGenerator; // generator of random jumps in BC
75 o2::math_utils::RandomRing<1000> mNCollBCGenerator; // generator of number of interactions in BC
76 o2::math_utils::RandomRing<1000> mCollTimeGenerator; // generator of number of interactions in BC
77
80 int mIntBCCache = 0;
81
82 float mIntRate = -1.;
83 float mBCTimeRMS = 0.2;
84 double mMuBC = -1.;
85
87 std::vector<float> mTimeInBC;
88 std::vector<uint16_t> mInteractingBCs; // vector of interacting BCs
89 int mCurrBCIdx = 0;
90
91 static constexpr float DefIntRate = 50e3;
92
94};
95
96//_________________________________________________
97inline void InteractionSampler::generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest)
98{
99 // fill vector with interaction records
100 dest.clear();
101 for (int i = dest.capacity(); i--;) {
102 dest.push_back(generateCollisionTime());
103 }
104}
105
106//_________________________________________________
108{
110 if ((mCurrBCIdx += n) >= (int)mInteractingBCs.size()) {
112 mCurrBCIdx %= mInteractingBCs.size();
113 }
115}
116
117// Special case of InteractionSampler without actual sampling.
118// Engineers interaction sequence by putting one in each N-th BC with multiplicity mult.
120{
121
122 public:
123 FixedSkipBC_InteractionSampler(int every_n, int mult) : mEveryN{every_n}, mMultiplicity{mult}, InteractionSampler() {}
124
125 protected:
126 int simulateInteractingBC() override;
127
128 private:
129 int mEveryN; // the skip number ---> fills every N-th BC in the bunch filling scheme
130 int mMultiplicity; // how many events to put if bc is filled
132};
133
134// A version of the interaction sampler which can sample according to non-uniform mu(bc) as
135// observed during data taking.
137{
138 public:
140 bool setBCIntensityScales(const std::vector<float>& scales_from_vector);
141 bool setBCIntensityScales(const TH1F& scales_from_histo); // initialize scales
142
143 // helper function to determine the scales from a histogram (count from event selection analysis)
144 std::vector<float> determineBCIntensityScalesFromHistogram(const TH1F& scales_from_histo);
145
146 const std::vector<float>& getBCIntensityScales() const { return mBCIntensityScales; }
147
148 protected:
149 int simulateInteractingBC() override;
150 int getBCJump() const;
151
152 private:
153 // non-uniformity
154 std::vector<float> mBCIntensityScales;
156};
157
158} // namespace steer
159} // namespace o2
160
161#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
ClassDef(InteractionSampler, 1)
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)
float mIntRate
total interaction rate in Hz
int mIntBCCache
N interactions left for current BC.
void setInteractionRate(float rateHz)
std::vector< float > determineBCIntensityScalesFromHistogram(const TH1F &scales_from_histo)
bool setBCIntensityScales(const std::vector< float > &scales_from_vector)
const std::vector< float > & getBCIntensityScales() const
GLdouble n
Definition glcorearb.h:1982
constexpr int LHCMaxBunches
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)