Belle II Software  release-05-02-19
SADInputModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010-2012 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Andreas Moll, Hiroyuki Nakayama *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <generators/modules/SADInputModule.h>
12 #include <mdst/dataobjects/MCParticleGraph.h>
13 
14 #include <framework/logging/Logger.h>
15 #include <framework/utilities/FileSystem.h>
16 
17 #include <framework/gearbox/Unit.h>
18 #include <framework/gearbox/GearDir.h>
19 
20 #include <framework/datastore/StoreArray.h>
21 
22 #include <TGeoMatrix.h>
23 
24 using namespace std;
25 using namespace Belle2;
26 
27 //-----------------------------------------------------------------
28 // Register the Module
29 //-----------------------------------------------------------------
30 REG_MODULE(SADInput)
31 
32 //-----------------------------------------------------------------
33 // Implementation
34 //-----------------------------------------------------------------
35 
37 {
38  //Set module properties
39  setDescription("Reads the SAD data from a root file and stores it into the MCParticle collection.");
40  setPropertyFlags(c_Input);
41 
42  //Parameter definition
43  addParam("AccRing", m_accRing, "The accelerator ring: 0 = LER, 1 = HER");
44  addParam("ReadoutTime", m_readoutTime, "The readout time of the detector [ns]", 20 * Unit::us);
45  addParam("ReadMode", m_readMode,
46  "The read mode: 0 = one SAD particle per event, 1 = one real particle per event, 2 = all SAD particles per event", 0);
47  addParam("Filename", m_filename, "The filename of the SAD input file.");
48  addParam("Range", m_range, "All particles within the range around the IP are loaded [cm].", 300.0 * Unit::cm);
49  addParam("PxResolution", m_pxRes, "The resolution for the x momentum component of the SAD real particle.", 0.01);
50  addParam("PyResolution", m_pyRes, "The resolution for the y momentum component of the SAD real particle.", 0.01);
51  addParam("RotateParticles", m_rotateParticles,
52  "Rotate the SAD particles around the nominal beam axis [deg] (just for unphysical tests !!!).", 0.0);
53 
54  m_PipePartMatrix = NULL;
55 }
56 
57 
58 void SADInputModule::initialize()
59 {
60  //Register inputs/outputs
61  m_eventMetaDataPtr.isRequired();
62  StoreArray<MCParticle> mcParticles;
63  mcParticles.registerInDataStore();
64 
65  //Check parameters
66  if (!FileSystem::fileExists(m_filename)) {
67  B2ERROR("Parameter <Filename>: Could not open the file. The filename " << m_filename << " does not exist !");
68  } else m_reader.open(m_filename);
69 
70  //Initialize the SAD reader.
71  //Set the transformation from local SAD plane to global geant4 space.
72  m_PipePartMatrix = new TGeoHMatrix("SADTrafo");
73  m_PipePartMatrix->RotateZ(m_rotateParticles);
74 
75  GearDir ler("/Detector/SuperKEKB/LER/");
76  GearDir her("/Detector/SuperKEKB/HER/");
77 
78  switch (m_accRing) {
79  case 0 : m_PipePartMatrix->RotateY(ler.getDouble("angle") / Unit::deg);
80  m_reader.initialize(m_PipePartMatrix, m_range, ReaderSAD::c_LER, m_readoutTime);
81  break;
82  case 1 : m_PipePartMatrix->RotateY(her.getDouble("angle") / Unit::deg);
83  m_reader.initialize(m_PipePartMatrix, m_range, ReaderSAD::c_HER, m_readoutTime);
84  break;
85  default: B2FATAL("Please specify a valid number for the accRing parameter (0 = LER, 1 = HER) !");
86  break;
87  }
88 
89  m_reader.setMomentumRes(m_pxRes, m_pyRes);
90 }
91 
92 
93 void SADInputModule::event()
94 {
95  try {
96  MCParticleGraph mpg;
97 
98  try {
99  //----------------------------------
100  // Read the SAD data
101  //----------------------------------
102  switch (m_readMode) {
103  case 0: readSADParticle(m_reader, mpg);
104  break;
105  case 1: readRealParticle(m_reader, mpg);
106  break;
107  case 2: m_reader.addAllSADParticles(mpg);
108  break;
109  default: B2FATAL("Please specify a valid number for the readMode parameter !");
110  break;
111  }
112 
113  //----------------------------------
114  // Generate MCParticles collection
115  //----------------------------------
116  mpg.generateList("", MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
117 
118  } catch (ReaderSAD::SADEndOfFile& exc) {
119  B2DEBUG(10, exc.what());
120  m_eventMetaDataPtr->setEndOfData();
121  return;
122  }
123  } catch (runtime_error& exc) {
124  B2ERROR(exc.what());
125  }
126 }
127 
128 
129 //====================================================================
130 // Private methods
131 //====================================================================
132 
133 void SADInputModule::readSADParticle(ReaderSAD& reader, MCParticleGraph& mpg)
134 {
135  double rate = reader.getSADParticle(mpg);
136  if (rate < 0) return;
137  m_eventMetaDataPtr->setGeneratedWeight(rate);
138 }
139 
140 void SADInputModule::readRealParticle(ReaderSAD& reader, MCParticleGraph& mpg)
141 {
142  if (!reader.getRealParticle(mpg)) return;
143  m_eventMetaDataPtr->setGeneratedWeight(1.0);
144 }
Belle2::MCParticleGraph::generateList
void generateList(const std::string &name="", int options=c_setNothing)
Generates the MCParticle list and stores it in the StoreArray with the given name.
Definition: MCParticleGraph.cc:211
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::MCParticleGraph
Class to build, validate and sort a particle decay chain.
Definition: MCParticleGraph.h:48
Belle2::SADInputModule
The SAD Input module.
Definition: SADInputModule.h:42
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::ProcType::c_Input
@ c_Input
Input Process.
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::GearDir
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:41
Belle2::ReaderSAD
Class to read files that have been created by SAD and store their content in a MCParticle graph.
Definition: ReaderSAD.h:45
Belle2::StoreArray< MCParticle >