Belle II Software  release-06-00-14
BeamBkgMixerModule.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <framework/core/Module.h>
12 #include <framework/datastore/StoreArray.h>
13 #include <simulation/background/BeamBGTypes.h>
14 #include <framework/dataobjects/BackgroundMetaData.h>
15 #include <string>
16 #include <map>
17 
18 #include "TChain.h"
19 #include "TClonesArray.h"
20 
21 
22 namespace Belle2 {
34  class BeamBkgMixerModule : public Module {
35 
36  public:
37 
42 
46  virtual ~BeamBkgMixerModule();
47 
52  virtual void initialize() override;
53 
58  virtual void beginRun() override;
59 
63  virtual void event() override;
64 
69  virtual void endRun() override;
70 
75  virtual void terminate() override;
76 
77  private:
78 
82  struct BkgHits {
83  TClonesArray* PXD;
84  TClonesArray* SVD;
85  TClonesArray* CDC;
86  TClonesArray* TOP;
87  TClonesArray* ARICH;
88  TClonesArray* ECL;
89  TClonesArray* BKLM;
90  TClonesArray* EKLM;
91  TClonesArray* BeamBackHits;
97  PXD(0), SVD(0), CDC(0), TOP(0), ARICH(0), ECL(0), BKLM(0), EKLM(0),
98  BeamBackHits(0)
99  {}
100  };
101 
105  struct BkgFiles {
107  std::string type;
108  double realTime;
109  double scaleFactor;
110  std::vector<std::string> fileNames;
112  std::unique_ptr<TChain> tree;
113  unsigned numFiles;
114  unsigned numEvents;
115  unsigned eventCount;
116  double rate;
117  unsigned index;
123  fileType(BackgroundMetaData::c_Usual),
124  tree(nullptr), numFiles(0), numEvents(0), eventCount(0), rate(0.0), index(0)
125  {}
137  const std::string& bkgType,
138  const std::string& fileName,
139  double time,
140  double scaleFac,
142  unsigned indx = 0):
143  tag(bkgTag), type(bkgType), realTime(time), scaleFactor(scaleFac),
144  fileType(fileTyp),
145  tree(nullptr), numFiles(0), numEvents(0), eventCount(0), rate(0.0), index(indx)
146  {
147  fileNames.push_back(fileName);
148  }
149  };
150 
151 
160  template<class SIMHIT>
162  TClonesArray* cloneArray,
163  double timeShift,
164  double minTime,
165  double maxTime)
166  {
167  if (!cloneArray) return;
168  if (!simHits.isValid()) return;
169 
170  int numEntries = cloneArray->GetEntriesFast();
171  for (int i = 0; i < numEntries; i++) {
172  SIMHIT* bkgSimHit = static_cast<SIMHIT*>(cloneArray->AddrAt(i));
173  SIMHIT* simHit = simHits.appendNew(*bkgSimHit);
174  simHit->shiftInTime(timeShift);
175  if (simHit->getBackgroundTag() == 0) // should be properly set at bkg simulation
176  simHit->setBackgroundTag(BackgroundMetaData::bg_other);
177  if (m_wrapAround) {
178  double time = simHit->getGlobalTime();
179  if (time > maxTime) {
180  double windowSize = maxTime - minTime;
181  double shift = int((time - minTime) / windowSize) * windowSize;
182  simHit->shiftInTime(-shift);
183  }
184  }
185  }
186 
187  }
188 
197  template<class HIT>
198  void addBeamBackHits(StoreArray<HIT>& hits, TClonesArray* cloneArray,
199  double timeShift, double minTime, double maxTime)
200  {
201  // Match SubDet id from BeamBackHits to whether we keep it or not
202  // KLM is a single component, but it has different hits for BKLM and EKLM.
203  bool keep[] = {false, m_PXD, m_SVD, m_CDC, m_ARICH, m_TOP, m_ECL, m_KLM, m_KLM};
204  if (!cloneArray) return;
205  if (!hits.isValid()) return;
206  // this is basically a copy of addSimHits but we only add the
207  // BeamBackHits from the specified sub detectors so we have to check
208  // each if it is from one of the enabled subdetectors
209  int numEntries = cloneArray->GetEntriesFast();
210  for (int i = 0; i < numEntries; i++) {
211  HIT* bkgHit = static_cast<HIT*>(cloneArray->AddrAt(i));
212  //Only keep selected
213  if (!keep[bkgHit->getSubDet()]) continue;
214  HIT* hit = hits.appendNew(*bkgHit);
215  hit->shiftInTime(timeShift);
216  //TODO: BeamBackHits does not have a setBackgroundTag so we do not
217  //check or set it
218  if (m_wrapAround) {
219  double time = hit->getTime();
220  if (time > maxTime) {
221  double windowSize = maxTime - minTime;
222  double shift = int((time - minTime) / windowSize) * windowSize;
223  hit->shiftInTime(-shift);
224  }
225  }
226  }
227  }
228 
236  bool isComponentIncluded(std::vector<std::string>& components,
237  const std::string& component);
238 
248  const std::string& bkgType,
249  const std::string& fileName,
250  double realTime,
252 
258  bool acceptEvent(TClonesArray* cloneArrayECL);
259 
260 
261  std::vector<std::string> m_backgroundFiles;
263  std::vector<std::tuple<std::string, double> > m_scaleFactors;
264  double m_minTime;
265  double m_maxTime;
266  std::vector<std::string> m_components;
268  double m_minTimeECL;
269  double m_maxTimeECL;
270  double m_minTimePXD;
271  double m_maxTimePXD;
272  double m_maxEdepECL;
275  std::vector<BkgFiles> m_backgrounds;
278  bool m_PXD = false;
279  bool m_SVD = false;
280  bool m_CDC = false;
281  bool m_TOP = false;
282  bool m_ARICH = false;
283  bool m_ECL = false;
284  bool m_KLM = false;
285  bool m_BeamBackHits = false;
289  std::map<std::string, int> m_rejected;
290  std::map<std::string, int> m_reused;
291  int m_rejectedCount = 0;
293  };
294 
296 } // Belle2 namespace
297 
Metadata information about the beam background file.
EFileType
Enum for BG file types.
BG_TAG
Enum for background tags.
@ bg_other
Other type of background.
New beam background mixer; this one doesn't need ROF files.
std::map< std::string, int > m_reused
messages: rejused events
double m_maxEdepECL
maximal allowed deposited energy in ECL
std::vector< BkgFiles > m_backgrounds
container for background samples
double m_maxTimeECL
maximal time shift of background event for ECL
std::map< std::string, int > m_rejected
messages: rejected events
bool m_ECL
true if found in m_components
std::vector< std::string > m_components
detector components
void addBeamBackHits(StoreArray< HIT > &hits, TClonesArray *cloneArray, double timeShift, double minTime, double maxTime)
functions that add BeamBackHits to those in the DataStore
double m_minTimePXD
minimal time shift of background event for PXD
bool m_PXD
true if found in m_components
void addSimHits(StoreArray< SIMHIT > &simHits, TClonesArray *cloneArray, double timeShift, double minTime, double maxTime)
functions that add background SimHits to those in the DataStore
background::BeamBGTypes m_bgTypes
defined BG types
std::vector< std::string > m_backgroundFiles
names of beam background files
bool m_CDC
true if found in m_components
double m_maxTime
maximal time shift of background event
double m_minTime
minimal time shift of background event
bool m_wrapAround
if true wrap around events in the tail after maxTime
double m_maxTimePXD
maximal time shift of background event for PXD
int m_rejectedCount
counter for suppresing "rejected event" messages
BkgHits m_simHits
input event buffer
int m_cacheSize
file cache size in Mbytes
bool m_BeamBackHits
if true add also background hits
bool m_ARICH
true if found in m_components
double m_overallScaleFactor
overall scale factor
double m_minTimeECL
minimal time shift of background event for ECL
std::vector< std::tuple< std::string, double > > m_scaleFactors
scale factors
bool m_KLM
true if found in m_components
bool m_SVD
true if found in m_components
bool m_TOP
true if found in m_components
Base class for Modules.
Definition: Module.h:72
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
Class to define BG types and to convert between BG types and tags or v.v.
Definition: BeamBGTypes.h:27
virtual ~BeamBkgMixerModule()
Destructor.
bool isComponentIncluded(std::vector< std::string > &components, const std::string &component)
Returns true if a component is found in components or the list is empty.
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
void appendSample(BackgroundMetaData::BG_TAG bkgTag, const std::string &bkgType, const std::string &fileName, double realTime, BackgroundMetaData::EFileType fileTyp)
appends background sample to m_backgrounds
bool acceptEvent(TClonesArray *cloneArrayECL)
Checks for deposited energy of ECLHits and returns true if Edep < m_maxEdepECL.
Abstract base class for different kinds of events.
structure to hold samples of a particular background type
std::vector< std::string > fileNames
file names
unsigned index
index of this element in the std::vector
double rate
background rate of the sample
std::unique_ptr< TChain > tree
tree pointer
BkgFiles(BackgroundMetaData::BG_TAG bkgTag, const std::string &bkgType, const std::string &fileName, double time, double scaleFac, BackgroundMetaData::EFileType fileTyp=BackgroundMetaData::c_Usual, unsigned indx=0)
usefull constructor
BackgroundMetaData::BG_TAG tag
background tag
double scaleFactor
scale factor for the rate
double realTime
real time of BG samlpe
unsigned numEvents
number of events (tree entries) in the sample
unsigned numFiles
number of files connected to TChain
BackgroundMetaData::EFileType fileType
file type
unsigned eventCount
current event (tree entry)
An input event buffer definition for background SimHits.
TClonesArray * BeamBackHits
BeamBackHits from collision file.
TClonesArray * ARICH
ARICH SimHits from collision file.
TClonesArray * BKLM
BKLM SimHits from collision file.
TClonesArray * CDC
CDC SimHits from collision file.
TClonesArray * PXD
PXD SimHits from collision file.
TClonesArray * ECL
ECL SimHits from collision file.
TClonesArray * SVD
SVD SimHits from collision file.
TClonesArray * EKLM
EKLM SimHits from collision file.
TClonesArray * TOP
TOP SimHits from collision file.