Belle II Software  release-08-01-10
ARICHReconstructorModule.cc
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 // Own header.
9 
10 #include <arich/modules/arichReconstruction/ARICHReconstructorModule.h>
11 #include <time.h>
12 
13 // Hit classes
14 #include <mdst/dataobjects/Track.h>
15 #include <tracking/dataobjects/ExtHit.h>
16 #include <mdst/dataobjects/MCParticle.h>
17 #include <arich/dataobjects/ARICHAeroHit.h>
18 #include <arich/dataobjects/ARICHLikelihood.h>
19 #include <arich/dataobjects/ARICHTrack.h>
20 #include <arich/dataobjects/ARICHHit.h>
21 
22 // framework - DataStore
23 #include <framework/datastore/DataStore.h>
24 #include <framework/datastore/StoreArray.h>
25 
26 // framework aux
27 #include <framework/gearbox/Unit.h>
28 #include <framework/gearbox/Const.h>
29 #include <framework/logging/Logger.h>
30 
31 // ROOT
32 #include <TVector3.h>
33 
34 using namespace boost;
35 
36 namespace Belle2 {
42  //-----------------------------------------------------------------
44  //-----------------------------------------------------------------
45 
46  REG_MODULE(ARICHReconstructor);
47 
48  //-----------------------------------------------------------------
49  // Implementation
50  //-----------------------------------------------------------------
51 
52  ARICHReconstructorModule::ARICHReconstructorModule() :
53  m_ana(NULL)
54  {
55  // Set description()
56  setDescription("This module calculates the ARICHLikelihood values for all particle id. hypotheses, for all tracks that enter ARICH in the event.");
57 
58  // Set property flags
60 
61  // Add parameters
62  addParam("trackPositionResolution", m_trackPositionResolution,
63  "Resolution of track position on aerogel plane (for additional smearing of MC tracks)", 1.0 * Unit::mm);
64  addParam("trackAngleResolution", m_trackAngleResolution,
65  "Resolution of track direction angle on aerogel plane (for additional smearing of MC tracks)", 2.0 * Unit::mrad);
66  addParam("inputTrackType", m_inputTrackType, "Input tracks switch: tracking (0), from AeroHits - MC info (1)", 0);
67  addParam("storePhotons", m_storePhot, "Set to 1 to store reconstructed photon information (Ch. angle,...)", 0);
68  addParam("useAlignment", m_align, "Use ARICH global position alignment constants", true);
69  addParam("useMirrorAlignment", m_alignMirrors, "Use ARICH mirror alignment constants", true);
70  }
71 
73  {
74  if (m_ana) delete m_ana;
75  }
76 
78  {
79  // Initialize variables
80  if (m_ana) delete m_ana;
85 
86  // Input: ARICHDigits
87  m_ARICHHits.isRequired();
88 
89  m_Tracks.isOptional();
90  m_ExtHits.isOptional();
91  m_aeroHits.isOptional();
92 
93  if (!m_aeroHits.isOptional()) {
94  m_Tracks.isRequired();
95  m_ExtHits.isRequired();
96  }
97 
98  StoreArray<MCParticle> MCParticles;
99  MCParticles.isOptional();
100 
101  // Output - log likelihoods
102  m_ARICHLikelihoods.registerInDataStore();
103 
104  m_ARICHTracks.registerInDataStore();
105 
106  if (m_inputTrackType) m_ARICHTracks.registerRelationTo(m_ARICHLikelihoods);
107  else {
108  m_ARICHTracks.registerRelationTo(m_ExtHits);
109  m_Tracks.registerRelationTo(m_ARICHLikelihoods);
110  }
111  //m_ARICHTracks.registerInDataStore(DataStore::c_DontWriteOut);
112  //m_ARICHTracks.registerRelationTo(m_ARICHLikelihoods, DataStore::c_Event, DataStore::c_DontWriteOut);
113  m_ARICHTracks.registerRelationTo(m_aeroHits);
115  }
116 
118  {
119  m_ana->initialize();
120  }
121 
123  {
124  // using track information form tracking system (mdst Track)
125  if (m_inputTrackType == 0) {
126 
127  Const::EDetector myDetID = Const::EDetector::ARICH; // arich
128  Const::ChargedStable hypothesis = Const::pion;
129  int pdgCode = abs(hypothesis.getPDGCode());
130 
131  for (const Track& track : m_Tracks) {
132  const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(hypothesis);
133  if (!fitResult) {
134  B2ERROR("No TrackFitResult for " << hypothesis.getPDGCode());
135  continue;
136  }
137 
138  const MCParticle* mcParticle = track.getRelated<MCParticle>();
139 
140  ARICHAeroHit* aeroHit = NULL;
141  if (mcParticle) aeroHit = mcParticle->getRelated<ARICHAeroHit>();
142 
143  RelationVector<ExtHit> extHits = DataStore::getRelationsWithObj<ExtHit>(&track);
144  ARICHTrack* arichTrack = NULL;
145 
146  //const ExtHit* arich1stHit = NULL;
147  const ExtHit* arich2ndHit = NULL;
148  const ExtHit* arichWinHit = NULL;
149 
150  for (unsigned i = 0; i < extHits.size(); i++) {
151  const ExtHit* extHit = extHits[i];
152  if (abs(extHit->getPdgCode()) != pdgCode or
153  extHit->getDetectorID() != myDetID or
154  extHit->getStatus() != EXT_EXIT or // particles registered at the EXIT of the Al plate
155  extHit->getMomentum().Z() < 0.0 or // track passes in backward
156  extHit->getCopyID() == 12345) { continue;}
157  if (extHit->getCopyID() == 6789) { arich2ndHit = extHit; continue;}
158  arichWinHit = extHit;
159  }
160 
161  if (arich2ndHit) {
162  // if aeroHit cannot be found using MCParticle check if it was already related to the extHit (by ARICHRelate module)
163  if (!aeroHit) aeroHit = arich2ndHit->getRelated<ARICHAeroHit>();
164 
165  // make new ARICHTrack
166  arichTrack = m_ARICHTracks.appendNew(arich2ndHit);
167  if (arichWinHit) arichTrack->setHapdWindowHit(arichWinHit);
168  }
169 
170  // skip if track has no extHit in ARICH
171  if (!arichTrack) continue;
172  // transform track parameters to ARICH local frame
173  m_ana->transformTrackToLocal(*arichTrack, m_align);
174  // make new ARICHLikelihood
175  ARICHLikelihood* like = m_ARICHLikelihoods.appendNew();
176  // calculate and set likelihood values
177  m_ana->likelihood2(*arichTrack, m_ARICHHits, *like);
178  // make relations
179  track.addRelationTo(like);
180  arichTrack->addRelationTo(arich2ndHit);
181  if (aeroHit) arichTrack->addRelationTo(aeroHit);
182 
183  } // Tracks loop
184  } // input type if
185 
186 
187  // using track information form MC (stored in ARICHAeroHit)
188  else {
189 
190  // Loop over all ARICHAeroHits
191  for (const ARICHAeroHit& aeroHit : m_aeroHits) {
192 
193  // make new ARICHTrack
194  ARICHTrack* arichTrack = m_ARICHTracks.appendNew(&aeroHit);
195  // smearing of track parameters (to mimic tracking system resolutions)
196  m_ana->smearTrack(*arichTrack);
197  // transform track parameters to ARICH local frame
198  m_ana->transformTrackToLocal(*arichTrack, m_align);
199  // make associated ARICHLikelihood
200  ARICHLikelihood* like = m_ARICHLikelihoods.appendNew();
201  // calculate and set likelihood values
202  m_ana->likelihood2(*arichTrack, m_ARICHHits, *like);
203  // make relation
204  arichTrack->addRelationTo(like);
205  arichTrack->addRelationTo(&aeroHit);
206  } // for iTrack
207 
208  }
209  }
210 
212  {
213  if (m_inputTrackType == 0) { B2DEBUG(100, "ARICHReconstructorModule: track infromation is taken from mdst Tracks.");}
214  else B2DEBUG(100, "ARICHReconstructorModule: track information is taken from MC (ARICHAeroHit).");
215  }
216 
218 } // namespace Belle2
Datastore class that holds information on track parameters at the entrance in aerogel.
Definition: ARICHAeroHit.h:27
This is a class to store ARICH likelihoods in the datastore.
Internal ARICH track reconstruction.
void useMirrorAlignment(bool align)
use mirror alignment or not
ARICHReconstruction * m_ana
Class with reconstruction tools.
bool m_alignMirrors
If==1 alignment constants are used for global->local track transformation.
int m_inputTrackType
Input tracks from the tracking (0) or from MCParticles>AeroHits (1).
double m_trackAngleResolution
Track direction resolution; simulation smearing.
int m_storePhot
If == 1 individual reconstruced photon information (cherenkov angle,...) is stored in ARICHTrack.
double m_trackPositionResolution
Track position resolution; simulation smearing.
Datastore class that holds position and momentum information of tracks that hit ARICH.
Definition: ARICHTrack.h:32
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:580
int getPDGCode() const
PDG code.
Definition: Const.h:464
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:32
int getPdgCode() const
Get PDG code of this extrapolation's hypothesis.
Definition: ExtHit.h:118
ExtHitStatus getStatus() const
Get state of extrapolation at this hit.
Definition: ExtHit.h:130
int getCopyID() const
Get detector-element ID of sensitive element within detector.
Definition: ExtHit.h:126
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
Definition: ExtHit.h:122
ROOT::Math::XYZVector getMomentum() const
Get momentum at this extrapolation hit.
Definition: ExtHit.h:159
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Class for type safe access to objects that are referred to in relations.
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
static const double mm
[millimeters]
Definition: Unit.h:70
static const double mrad
[millirad]
Definition: Unit.h:108
void setTrackPositionResolution(double pRes)
Sets track position resolution (from tracking)
int likelihood2(ARICHTrack &arichTrack, const StoreArray< ARICHHit > &arichHits, ARICHLikelihood &arichLikelihood)
Computes the value of identity likelihood function for different particle hypotheses.
void initialize()
read geomerty parameters from xml and initialize class memebers
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
void setTrackAngleResolution(double aRes)
Sets track direction resolution (from tracking)
virtual void beginRun() override
Called when entering a new run.
void printModuleParams()
Print module parameters.
void transformTrackToLocal(ARICHTrack &arichTrack, bool align)
Transforms track parameters from global Belle2 to ARICH local frame.
int smearTrack(ARICHTrack &arichTrack)
Smeares track parameters ("simulate" the uncertainties of tracking).
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.