10#include <ecl/modules/eclHadronTimeCalibrationValidationCollector/eclHadronTimeCalibrationValidationCollectorModule.h>
13#include <ecl/dataobjects/ECLCalDigit.h>
14#include <ecl/dataobjects/ECLElementNumbers.h>
17#include <framework/gearbox/Const.h>
18#include <mdst/dataobjects/ECLCluster.h>
19#include <mdst/dataobjects/Track.h>
20#include <mdst/dataobjects/HitPatternCDC.h>
32REG_MODULE(eclHadronTimeCalibrationValidationCollector);
40 m_dbg_tree_photonClusters(0),
48 "Events with fabs(getTimeFit) > m_timeAbsMax "
49 "are excluded", (
short)80);
52 "If true, TTree 'tree' with more detailed event info is saved in "
53 "the output file specified by HistoManager",
77 "Validating crystal and crate time calibrations using photon clusters in events with lots of tracks and clusters") ;
93 "Validating crystal and crate time calibrations using photon clusters in events with lots of tracks and clusters") ;
110 B2INFO(
"eclHadronTimeCalibrationValidationCollector: Experiment = " <<
m_EventMetaData->getExperiment() <<
116 double binSize = 2000.0 / pow(2, 12);
117 double halfBinSize = binSize / 2.0;
121 double nBinsNeg = floor((
m_timeAbsMax - halfBinSize) / binSize);
122 double min_t = -nBinsNeg * binSize - halfBinSize;
123 int nbins = nBinsNeg * 2 + 1;
124 double max_t = min_t + nbins * binSize;
128 const Int_t N_E_BIN_EDGES = 64;
129 const Int_t N_E_BINS = N_E_BIN_EDGES - 1;
130 Double_t energyBinEdges[N_E_BIN_EDGES] = {0, 0.05, 0.051, 0.052, 0.053, 0.054, 0.055, 0.056, 0.057, 0.058, 0.059, 0.06, 0.062, 0.064, 0.066, 0.068, 0.07, 0.075, 0.08, 0.085, 0.09, 0.095, 0.1, 0.11, 0.12, 0.13, 0.14, 0.15, 0.16, 0.17, 0.18, 0.19, 0.2, 0.25, 0.3, 0.35, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1, 1.2, 1.4, 1.6, 1.8, 2, 2.25, 2.5, 2.8, 3.2, 3.6, 4, 4.4, 4.8, 5.2, 5.6, 6, 6.4, 6.8, 7.2, 7.6, 8};
133 auto cutflow =
new TH1F(
"cutflow",
" ;Cut label number ;Number of events passing cut", 10, 0, 10) ;
134 registerObject<TH1F>(
"cutflow", cutflow) ;
136 auto clusterTime =
new TH1F(
"clusterTime",
";Photon ECL cluster time [ns]; number of photon ECL clusters", nbins, min_t, max_t) ;
137 registerObject<TH1F>(
"clusterTime", clusterTime) ;
139 auto clusterTime_cid =
new TH2F(
"clusterTime_cid",
142 registerObject<TH2F>(
"clusterTime_cid", clusterTime_cid) ;
144 auto clusterTime_run =
new TH2F(
"clusterTime_run",
145 ";Run number ;Photon ECL cluster time [ns]", 7000, 0, 7000, nbins, min_t, max_t) ;
146 registerObject<TH2F>(
"clusterTime_run", clusterTime_run) ;
149 auto clusterTimeClusterE =
new TH2F(
"clusterTimeClusterE",
150 ";Photon cluster energy [GeV];Photon cluster time [ns]", N_E_BINS, energyBinEdges, nbins, min_t, max_t) ;
151 registerObject<TH2F>(
"clusterTimeClusterE", clusterTimeClusterE) ;
154 auto dt99_clusterE =
new TH2F(
"dt99_clusterE",
155 ";Photon cluster energy [GeV];dt99 [ns]", N_E_BINS, energyBinEdges, nbins, 0, max_t) ;
156 registerObject<TH2F>(
"dt99_clusterE", dt99_clusterE) ;
159 auto eventT0 =
new TH1F(
"eventT0",
";event t0 [ns]; number of events", nbins, min_t, max_t) ;
160 registerObject<TH1F>(
"eventT0", eventT0) ;
163 auto clusterTimeE0E1diff =
new TH1F(
"clusterTimeE0E1diff",
164 ";ECL cluster time of max E photon - ECL cluster time of 2nd max E photon [ns]; number of photon ECL cluster time differences",
165 nbins, min_t, max_t) ;
166 registerObject<TH1F>(
"clusterTimeE0E1diff", clusterTimeE0E1diff) ;
178 int cutIndexPassed = 0;
179 getObjectPtr<TH1F>(
"cutflow")->Fill(cutIndexPassed);
180 B2DEBUG(22,
"Cutflow: no cuts: index = " << cutIndexPassed);
187 int tempCrysID = eclCalDigit.getCellId() - 1;
188 m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
193 double evt_t0 = -1000 ;
194 double evt_t0_unc = -1000 ;
206 evt_t0_unc =
m_eventT0->getEventT0Uncertainty() ;
208 B2DEBUG(26,
"Found event t0") ;
214 int nTrkAll =
tracks.getEntries() ;
223 for (
int iTrk = 0 ; iTrk < nTrkAll ; iTrk++) {
227 if (not tempTrackFit) {continue ;}
231 double z0 = tempTrackFit->
getZ0() ;
232 double d0 = tempTrackFit->
getD0() ;
271 B2DEBUG(26,
"Found loose and tight tracks") ;
274 int numGoodTightTracks_minCut = 4 ;
275 if (nTrkTight < numGoodTightTracks_minCut) {
280 getObjectPtr<TH1F>(
"cutflow")->Fill(cutIndexPassed) ;
281 B2DEBUG(22,
"Cutflow: At least " << numGoodTightTracks_minCut <<
" tight tracks: index = " << cutIndexPassed) ;
284 int numGoodLooseTracks_minCut = numGoodTightTracks_minCut ;
285 if (nTrkLoose < numGoodLooseTracks_minCut) {
290 getObjectPtr<TH1F>(
"cutflow")->Fill(cutIndexPassed) ;
291 B2DEBUG(22,
"Cutflow: No additional loose tracks: index = " << cutIndexPassed) ;
296 double clusterE_minCut = 0.1 ;
298 int nGoodClusts = 0 ;
299 vector<int> goodClusterIdxs ;
300 for (
int ic = 0; ic < nclust; ic++) {
302 if (eClust > clusterE_minCut) {
303 goodClusterIdxs.push_back(ic) ;
310 int numGoodEMclusters_minCut = 5 ;
311 if (nGoodClusts < numGoodEMclusters_minCut) {
316 getObjectPtr<TH1F>(
"cutflow")->Fill(cutIndexPassed) ;
317 B2DEBUG(22,
"Cutflow: At least " << numGoodEMclusters_minCut <<
" ECL clusters: index = " << cutIndexPassed) ;
323 double photonE_minCut = 0.05 ;
324 double zernikeMVA_minCut = 0.2 ;
327 vector<int> goodPhotonClusterIdxs ;
328 for (
int ic = 0; ic < nclust; ic++) {
336 if ((eClust > photonE_minCut) &&
339 (!badPhotonTimeResolution) &&
340 (zernikeMVA > zernikeMVA_minCut) &&
342 goodPhotonClusterIdxs.push_back(ic) ;
350 int numGoodPhotonclusters_minCut = 1 ;
351 if (nPhotons < numGoodPhotonclusters_minCut) {
356 getObjectPtr<TH1F>(
"cutflow")->Fill(cutIndexPassed) ;
357 B2DEBUG(22,
"Cutflow: At least " << numGoodPhotonclusters_minCut <<
" good photon: index = " << cutIndexPassed) ;
364 vector<double> goodClustTimes ;
365 vector<double> goodClust_dt99 ;
366 vector<double> goodClusters_crysE ;
367 vector<double> goodClustE ;
368 vector<int> goodClustMaxEcrys_cid ;
369 for (
long unsigned int i = 0; i < goodPhotonClusterIdxs.size(); i++) {
370 int ic = goodPhotonClusterIdxs[i] ;
376 goodClustMaxEcrys_cid.push_back(cid) ;
379 goodClusters_crysE.push_back(
m_EperCrys[cid - 1]) ;
380 goodClustE.push_back(eClust);
387 vector< pair<double, double> > pair_energy_time ;
388 for (
long unsigned int ic = 0; ic < goodClusters_crysE.size(); ic++) {
389 pair_energy_time.push_back(make_pair(goodClusters_crysE[ic], goodClustTimes[ic])) ;
394 sort(pair_energy_time.begin(), pair_energy_time.end(), greater<>()) ;
398 B2DEBUG(22,
"Event passed all cuts");
402 getObjectPtr<TH1F>(
"eventT0")->Fill(evt_t0) ;
404 bool isCDCt0 =
m_eventT0->isCDCEventT0();
405 bool isECLt0 =
m_eventT0->isECLEventT0();
406 string t0Detector =
"UNKNOWN... WHY?";
409 }
else if (isECLt0) {
413 B2DEBUG(26,
"t0 = " << evt_t0 <<
" ns. t0 is from CDC?=" << isCDCt0 <<
", t0 is from ECL?=" << isECLt0 <<
" t0 from " <<
419 for (
long unsigned int i = 0 ; i < goodPhotonClusterIdxs.size() ; i++) {
420 getObjectPtr<TH1F>(
"clusterTime")->Fill(goodClustTimes[i]) ;
421 getObjectPtr<TH2F>(
"clusterTime_cid")->Fill(goodClustMaxEcrys_cid[i] + 0.001, goodClustTimes[i], 1) ;
422 getObjectPtr<TH2F>(
"clusterTime_run")->Fill(
m_EventMetaData->getRun() + 0.001, goodClustTimes[i], 1) ;
423 getObjectPtr<TH2F>(
"clusterTimeClusterE")->Fill(goodClustE[i], goodClustTimes[i], 1) ;
424 getObjectPtr<TH2F>(
"dt99_clusterE")->Fill(goodClustE[i], goodClust_dt99[i], 1) ;
445 B2DEBUG(26,
"Filled cluster tree") ;
448 if (pair_energy_time.size() >= 2) {
449 getObjectPtr<TH1F>(
"clusterTimeE0E1diff")->Fill(pair_energy_time[0].second - pair_energy_time[1].second) ;
468 B2DEBUG(26,
"Filled event tree") ;
Calibration collector module base class.
static const ChargedStable pion
charged pion particle
@ c_nPhotons
CR is split into n photons (N1)
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Values of the result of a track fit with a given particle hypothesis.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
double m_E_photon_clust
Photon cluster energy.
StoreObjPtr< EventT0 > m_eventT0
StoreObjPtr for T0.
int m_tree_run
Run number for debug TTree output.
double m_looseTrkD0
Loose track d0 minimum cut.
bool m_saveTree
If true, save TTree with more detailed event info.
StoreArray< ECLCluster > m_eclClusterArray
Required input array of ECLClusters.
short m_timeAbsMax
Events with abs(time) > m_timeAbsMax are excluded, mostly for histogram x-range purposes.
int m_tree_cid
ECL Cell ID (1..ECLElementNumbers::c_NCrystals) for debug TTree output.
TTree * m_dbg_tree_event
debug output tree for per event
double m_tree_time
Calibrated time.
double m_tightTrkZ0
Tight track z0 minimum cut.
double m_tree_E0
Highest E cluster energy.
int m_NtightTracks
Number of tight tracks.
TTree * m_dbg_tree_photonClusters
Output tree with detailed event data.
void collect() override
Select events and crystals and accumulate histograms.
virtual ~eclHadronTimeCalibrationValidationCollectorModule()
Module destructor.
double m_tree_dt99
dt99 for cluster
StoreArray< ECLCalDigit > m_eclCalDigitArray
Required input array of ECLCalDigits.
std::vector< float > m_EperCrys
ECL Cal digit energy for each crystal.
void prepare() override
Define histograms and read payloads from DB.
double m_tightTrkD0
Tight track d0 minimum cut.
double m_tree_t0
EventT0 (not from ECL) for debug TTree output.
double m_tree_time_fromE0
Calibrated time - highest E cluster.
void inDefineHisto() override
Replacement for defineHisto() in CalibrationCollector modules.
int m_NphotonClusters
Number of photon clusters.
StoreArray< Track > tracks
Required input array of tracks.
eclHadronTimeCalibrationValidationCollectorModule()
Module constructor.
int m_NGoodClusters
Number of good clusters.
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
double m_tree_t0_unc
EventT0 uncertainty for debug TTree output.
int m_tree_evt_num
Event number for debug TTree output.
double m_looseTrkZ0
Loose track z0 minimum cut.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.