forked from AliceO2Group/AliceO2
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathInteractionSampler.h
More file actions
136 lines (113 loc) · 4.49 KB
/
InteractionSampler.h
File metadata and controls
136 lines (113 loc) · 4.49 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
/// @brief Simulated interaction record sampler
#ifndef ALICEO2_INTERACTIONSAMPLER_H
#define ALICEO2_INTERACTIONSAMPLER_H
#include <Rtypes.h>
#include <TMath.h>
#include <TRandom.h>
#include <vector>
#include "CommonDataFormat/InteractionRecord.h"
#include "CommonDataFormat/BunchFilling.h"
#include "CommonConstants/LHCConstants.h"
#include "MathUtils/RandomRing.h"
namespace o2
{
namespace steer
{
class InteractionSampler
{
public:
static constexpr float Sec2NanoSec = 1.e9; // s->ns conversion
const o2::InteractionTimeRecord& generateCollisionTime();
void generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest);
void init();
void setInteractionRate(float rateHz)
{
mIntRate = rateHz;
mMuBC = -1.; // invalidate
}
float getInteractionRate() const { return mIntRate; }
void setFirstIR(const o2::InteractionRecord& ir)
{
mFirstIR.InteractionRecord::operator=(ir);
if (mFirstIR.orbit == 0 && mFirstIR.bc < 4) {
mFirstIR.bc = 4;
}
}
const o2::InteractionRecord& getFirstIR() const { return mFirstIR; }
void setMuPerBC(float mu)
{
mMuBC = mu;
mIntRate = -1.; // invalidate
}
float getMuPerBC() const { return mMuBC; }
void setBCTimeRMS(float tNS = 0.2) { mBCTimeRMS = tNS; }
float getBCTimeRMS() const { return mBCTimeRMS; }
const BunchFilling& getBunchFilling() const { return mBCFilling; }
BunchFilling& getBunchFilling() { return mBCFilling; }
void setBunchFilling(const BunchFilling& bc) { mBCFilling = bc; }
void setBunchFilling(const std::string& bcFillingFile);
void print() const;
protected:
virtual int simulateInteractingBC();
void nextCollidingBC(int n);
o2::math_utils::RandomRing<10000> mBCJumpGenerator; // generator of random jumps in BC
o2::math_utils::RandomRing<1000> mNCollBCGenerator; // generator of number of interactions in BC
o2::math_utils::RandomRing<1000> mCollTimeGenerator; // generator of number of interactions in BC
o2::InteractionTimeRecord mIR{{0, 0}, 0.};
o2::InteractionTimeRecord mFirstIR{{4, 0}, 0.};
int mIntBCCache = 0; ///< N interactions left for current BC
float mIntRate = -1.; ///< total interaction rate in Hz
float mBCTimeRMS = 0.2; ///< BC time spread in NANOSECONDS
double mMuBC = -1.; ///< interaction probability per BC
o2::BunchFilling mBCFilling; ///< patter of active BCs
std::vector<float> mTimeInBC; ///< interaction times within single BC
std::vector<uint16_t> mInteractingBCs; // vector of interacting BCs
int mCurrBCIdx = 0; ///< counter for current interacting bunch
static constexpr float DefIntRate = 50e3; ///< default interaction rate
ClassDef(InteractionSampler, 1);
};
//_________________________________________________
inline void InteractionSampler::generateCollisionTimes(std::vector<o2::InteractionTimeRecord>& dest)
{
// fill vector with interaction records
dest.clear();
for (int i = dest.capacity(); i--;) {
dest.push_back(generateCollisionTime());
}
}
//_________________________________________________
inline void InteractionSampler::nextCollidingBC(int n)
{
/// get colliding BC as n-th after current one
if ((mCurrBCIdx += n) >= (int)mInteractingBCs.size()) {
mIR.orbit += mCurrBCIdx / mInteractingBCs.size();
mCurrBCIdx %= mInteractingBCs.size();
}
mIR.bc = mInteractingBCs[mCurrBCIdx];
}
// Special case of InteractionSampler without actual sampling.
// Engineers interaction sequence by putting one in each N-th BC with multiplicity mult.
class FixedSkipBC_InteractionSampler : public InteractionSampler
{
public:
FixedSkipBC_InteractionSampler(int every_n, int mult) : mEveryN{every_n}, mMultiplicity{mult}, InteractionSampler() {}
protected:
int simulateInteractingBC() override;
private:
int mEveryN; // the skip number ---> fills every N-th BC in the bunch filling scheme
int mMultiplicity; // how many events to put if bc is filled
ClassDef(FixedSkipBC_InteractionSampler, 1);
};
} // namespace steer
} // namespace o2
#endif