9 #include <ecl/modules/eclee5x5Collector/eclee5x5CollectorModule.h>
15 #include <analysis/ClusterUtility/ClusterUtils.h>
18 #include <framework/dataobjects/EventMetaData.h>
19 #include <framework/datastore/RelationVector.h>
22 #include <mdst/dataobjects/TRGSummary.h>
23 #include <mdst/dataobjects/ECLCluster.h>
26 #include <ecl/dataobjects/ECLDigit.h>
27 #include <ecl/dataobjects/ECLCalDigit.h>
28 #include <ecl/dbobjects/ECLCrystalCalib.h>
47 m_ECLExpee5x5E("ECLExpee5x5E"), m_ElectronicsCalib("ECLCrystalElectronics"), m_ee5x5Calib("ECLCrystalEnergyee5x5"),
48 m_selectdPhiData("ECLeedPhiData"), m_selectdPhiMC("ECLeedPhiMC")
51 setDescription(
"Calibration Collector Module for ECL single crystal energy calibration using Bhabha events");
52 setPropertyFlags(c_ParallelProcessingCertified);
53 addParam(
"thetaLabMinDeg", m_thetaLabMinDeg,
"miniumum ecl cluster theta in lab (degrees)", 17.);
54 addParam(
"thetaLabMaxDeg", m_thetaLabMaxDeg,
"maximum ecl cluster theta in lab (degrees)", 150.);
55 addParam(
"minE0", m_minE0,
"minimum energy of cluster 0: E*0/sqrts", 0.45);
56 addParam(
"minE1", m_minE1,
"minimum energy of cluster 1: E*1/sqrts", 0.40);
57 addParam(
"maxdThetaSum", m_maxdThetaSum,
"maximum diff between 180 deg and sum of cluster theta* (deg)", 2.);
58 addParam(
"dPhiScale", m_dPhiScale,
"scale dPhi* cut by this factor", 1.);
59 addParam(
"maxTime", m_maxTime,
"maximum cluster time diff abs(t1-t0)/dt99", 10.);
60 addParam(
"useCalDigits", m_useCalDigits,
"use MC events to obtain expected energies",
false);
61 addParam(
"requireL1", m_requireL1,
"only use events that have a level 1 trigger",
false);
68 void eclee5x5CollectorModule::prepare()
72 m_sqrts = m_boostrotate.getCMSEnergy();
73 B2INFO(
"eclee5x5Collector: Experiment = " << m_evtMetaData->getExperiment() <<
" run = " << m_evtMetaData->getRun() <<
" sqrts = "
79 auto EnVsCrysID =
new TH2F(
"EnVsCrysID",
"Normalized 5x5 energy for each crystal;crystal ID;E25/Expected", 8736, 0, 8736, 200, 0.7,
81 registerObject<TH2F>(
"EnVsCrysID", EnVsCrysID);
83 auto RvsCrysID =
new TH1F(
"RvsCrysID",
"E_exp x E_crysID / sigma^2;crysID;sum of E_exp x E_crysID/sigma^2", 8736, 0, 8736);
84 registerObject(
"RvsCrysID", RvsCrysID);
86 auto NRvsCrysID =
new TH1F(
"NRvsCrysID",
"Entries in RvsCrysID vs crysID;crysID;Entries in RvsCrysID", 8736, 0, 8736);
87 registerObject(
"NRvsCrysID", NRvsCrysID);
89 auto Qmatrix =
new TH2F(
"Qmatrix",
"E_i x E_j/sigma^2;crysID i;crysID j", 8736, 0, 8736, 8736, 0, 8736);
90 registerObject(
"Qmatrix", Qmatrix);
93 auto ElecCalibvsCrys =
new TH1F(
"ElecCalibvsCrys",
"Sum electronics calib const vs crystal ID;crystal ID;calibration constant",
95 registerObject<TH1F>(
"ElecCalibvsCrys", ElecCalibvsCrys);
96 auto ExpEvsCrys =
new TH1F(
"ExpEvsCrys",
"Sum expected energy calib const vs crystalID;crystal ID;calibration constant", 8736, 0,
98 registerObject<TH1F>(
"ExpEvsCrys", ExpEvsCrys);
99 auto InitialCalibvsCrys =
new TH1F(
"InitialCalibvsCrys",
"Sum initial calib const vs crystal ID;crystal ID;calibration constant",
101 registerObject<TH1F>(
"InitialCalibvsCrys", InitialCalibvsCrys);
103 auto CalibEntriesvsCrys =
new TH1F(
"CalibEntriesvsCrys",
"Entries in calib vs crys histograms;crystal ID;Entries per crystal", 8736,
105 registerObject<TH1F>(
"CalibEntriesvsCrys", CalibEntriesvsCrys);
107 auto EntriesvsCrys =
new TH1F(
"EntriesvsCrys",
"Selected Bhabha clusters vs crystal ID;crystal ID;Entries", 8736, 0, 8736);
108 registerObject<TH1F>(
"EntriesvsCrys", EntriesvsCrys);
110 auto dPhivsThetaID =
new TH2F(
"dPhivsThetaID",
111 "Phi* vs thetaID forward, pass thetaSum,E0,E1;thetaID of forward cluster;dPhi COM (deg)", 69, 0, 69, 150, 165, 180);
112 registerObject<TH2F>(
"dPhivsThetaID", dPhivsThetaID);
116 B2INFO(
"Input parameters to eclee5x5Collector:");
117 B2INFO(
"thetaLabMinDeg: " << m_thetaLabMinDeg);
118 B2INFO(
"thetaLabMaxDeg: " << m_thetaLabMaxDeg);
119 m_thetaLabMin = m_thetaLabMinDeg / TMath::RadToDeg();
120 m_thetaLabMax = m_thetaLabMaxDeg / TMath::RadToDeg();
121 B2INFO(
"minE0: " << m_minE0);
122 B2INFO(
"minE1: " << m_minE1);
123 B2INFO(
"maxdThetaSum: " << m_maxdThetaSum);
124 B2INFO(
"dPhiScale: " << m_dPhiScale);
125 B2INFO(
"maxTime: " << m_maxTime);
126 B2INFO(
"useCalDigits: " << m_useCalDigits);
127 B2INFO(
"requireL1: " << m_requireL1);
130 EperCrys.resize(8736);
131 m_thetaID.resize(8736);
138 if (m_ECLExpee5x5E.hasChanged()) {
139 Expee5x5E = m_ECLExpee5x5E->getCalibVector();
140 Expee5x5Sigma = m_ECLExpee5x5E->getCalibUncVector();
142 if (m_ElectronicsCalib.hasChanged()) {ElectronicsCalib = m_ElectronicsCalib->getCalibVector();}
143 if (m_ee5x5Calib.hasChanged()) {ee5x5Calib = m_ee5x5Calib->getCalibVector();}
144 if (m_selectdPhiMC.hasChanged() and m_useCalDigits) {
145 meandPhi = m_selectdPhiMC->getCalibVector();
146 widthdPhi = m_selectdPhiMC->getCalibUncVector();
147 }
else if (m_selectdPhiData.hasChanged()) {
148 meandPhi = m_selectdPhiData->getCalibVector();
149 widthdPhi = m_selectdPhiData->getCalibUncVector();
153 for (
int ic = 1; ic < 9000; ic += 1000) {
154 B2INFO(
"DB constants for cellID=" << ic <<
": ee5x5Calib = " << ee5x5Calib[ic - 1] <<
" Expee5x5E = " << Expee5x5E[ic - 1] <<
155 " ElectronicsCalib = " <<
156 ElectronicsCalib[ic - 1]);
160 for (
int crysID = 0; crysID < 8736; crysID++) {
161 if (ElectronicsCalib[crysID] <= 0) {B2FATAL(
"eclee5x5Collector: ElectronicsCalib = " << ElectronicsCalib[crysID] <<
" for crysID = " << crysID);}
162 if (Expee5x5E[crysID] == 0) {B2FATAL(
"eclee5x5Collector: Expee5x5E = 0 for crysID = " << crysID);}
163 if (ee5x5Calib[crysID] == 0) {B2FATAL(
"eclee5x5Collector: ee5x5Calib = 0 for crysID = " << crysID);}
168 m_eclClusterArray.isRequired();
169 m_eclCalDigitArray.isRequired();
170 m_eclDigitArray.isRequired();
171 m_evtMetaData.isRequired();
176 for (
int it = 0; it < 69; it++) {
179 m_dPhiMin.push_back(meandPhi[crysID] - m_dPhiScale * widthdPhi[crysID]);
180 m_dPhiMax.push_back(meandPhi[crysID] + m_dPhiScale * widthdPhi[crysID]);
181 for (
int ic = 0; ic < m_eclNeighbours5x5->getCrystalsPerRing(it); ic++) {
182 m_thetaID.at(crysID) = it;
192 void eclee5x5CollectorModule::collect()
197 for (
int crysID = 0; crysID < 8736; crysID++) {
198 getObjectPtr<TH1F>(
"ExpEvsCrys")->Fill(crysID + 0.001, Expee5x5E[crysID]);
199 getObjectPtr<TH1F>(
"ElecCalibvsCrys")->Fill(crysID + 0.001, ElectronicsCalib[crysID]);
200 getObjectPtr<TH1F>(
"InitialCalibvsCrys")->Fill(crysID + 0.001, ee5x5Calib[crysID]);
201 getObjectPtr<TH1F>(
"CalibEntriesvsCrys")->Fill(crysID + 0.001);
208 bool newConst =
false;
209 if (m_ECLExpee5x5E.hasChanged()) {
211 B2INFO(
"ECLExpee5x5E has changed, exp = " << m_evtMetaData->getExperiment() <<
" run = " << m_evtMetaData->getRun());
212 Expee5x5E = m_ECLExpee5x5E->getCalibVector();
213 Expee5x5Sigma = m_ECLExpee5x5E->getCalibUncVector();
215 if (m_ElectronicsCalib.hasChanged()) {
217 B2INFO(
"ECLCrystalElectronics has changed, exp = " << m_evtMetaData->getExperiment() <<
" run = " << m_evtMetaData->getRun());
218 ElectronicsCalib = m_ElectronicsCalib->getCalibVector();
220 if (m_ee5x5Calib.hasChanged()) {
222 B2INFO(
"ECLCrystalEnergyee5x5 has changed, exp = " << m_evtMetaData->getExperiment() <<
" run = " << m_evtMetaData->getRun());
223 ee5x5Calib = m_ee5x5Calib->getCalibVector();
227 for (
int ic = 1; ic < 9000; ic += 1000) {
228 B2INFO(
"DB constants for cellID=" << ic <<
": ee5x5Calib = " << ee5x5Calib[ic - 1] <<
" Expee5x5E = " << Expee5x5E[ic - 1] <<
229 " ElectronicsCalib = " <<
230 ElectronicsCalib[ic - 1]);
234 for (
int crysID = 0; crysID < 8736; crysID++) {
235 if (ElectronicsCalib[crysID] <= 0) {B2FATAL(
"eclee5x5Collector: ElectronicsCalib = " << ElectronicsCalib[crysID] <<
" for crysID = " << crysID);}
236 if (Expee5x5E[crysID] == 0) {B2FATAL(
"eclee5x5Collector: Expee5x5E = 0 for crysID = " << crysID);}
237 if (ee5x5Calib[crysID] == 0) {B2FATAL(
"eclee5x5Collector: ee5x5Calib = 0 for crysID = " << crysID);}
244 unsigned int L1TriggerResults = m_TRGResults->getTRGSummary(0);
245 if (L1TriggerResults == 0) {
return;}
250 int icMax[2] = { -1, -1};
251 double maxClustE[2] = { -1., -1.};
252 int nclust = m_eclClusterArray.getEntries();
253 for (
int ic = 0; ic < nclust; ic++) {
256 if (eClust > maxClustE[0]) {
257 maxClustE[1] = maxClustE[0];
259 maxClustE[0] = eClust;
261 }
else if (eClust > maxClustE[1]) {
262 maxClustE[1] = eClust;
271 if (icMax[1] == -1) {
return;}
272 double theta0 = m_eclClusterArray[icMax[0]]->getTheta();
273 double theta1 = m_eclClusterArray[icMax[1]]->getTheta();
274 if (theta0 < m_thetaLabMin || theta0 > m_thetaLabMax || theta1 < m_thetaLabMin || theta1 > m_thetaLabMax) {
return;}
277 double t0 = m_eclClusterArray[icMax[0]]->getTime();
278 double dt990 = m_eclClusterArray[icMax[0]]->getDeltaTime99();
279 double t1 = m_eclClusterArray[icMax[1]]->getTime();
280 double dt991 = m_eclClusterArray[icMax[1]]->getDeltaTime99();
281 double dt99min = dt990;
282 if (dt991 < dt990) {dt99min = dt991;}
283 if (dt99min <= 0) {dt99min = 0.0001;}
284 if (abs(t1 - t0) > dt99min * m_maxTime) {
return;}
291 double phi0 = m_eclClusterArray[icMax[0]]->getPhi();
292 TVector3 p30(0., 0., maxClustE[0]);
293 p30.SetTheta(theta0);
296 ECLCluster::EHypothesisBit::c_nPhotons);
298 double phi1 = m_eclClusterArray[icMax[1]]->getPhi();
299 TVector3 p31(0., 0., maxClustE[1]);
300 p31.SetTheta(theta1);
303 ECLCluster::EHypothesisBit::c_nPhotons);
306 TLorentzVector p40COM = m_boostrotate.rotateLabToCms() * p40;
307 TLorentzVector p41COM = m_boostrotate.rotateLabToCms() * p41;
308 double theta01COM = (p41COM.Theta() + p40COM.Theta()) * TMath::RadToDeg();
309 if (abs(theta01COM - 180.) > m_maxdThetaSum) {
return;}
312 if (p40COM.E() < m_minE0 * m_sqrts and p41COM.E() < m_minE0 * m_sqrts) {
return;}
313 if (p40COM.E() < m_minE1 * m_sqrts or p41COM.E() < m_minE1 * m_sqrts) {
return;}
318 int crysIDMax[2] = { -1, -1};
319 double crysEMax[2] = { -1., -1.};
320 for (
int imax = 0; imax < 2; imax++) {
321 auto eclClusterRelations = m_eclClusterArray[icMax[imax]]->getRelationsTo<
ECLCalDigit>(
"ECLCalDigits");
322 for (
unsigned int ir = 0; ir < eclClusterRelations.size(); ir++) {
323 const auto calDigit = eclClusterRelations.object(ir);
324 int tempCrysID = calDigit->
getCellId() - 1;
325 float tempE = calDigit->getEnergy();
326 if (tempE > crysEMax[imax]) {
327 crysEMax[imax] = tempE;
328 crysIDMax[imax] = tempCrysID;
334 double dphiCOM = abs(p41COM.Phi() - p40COM.Phi()) * TMath::RadToDeg();
335 if (dphiCOM > 180.) {dphiCOM = 360. - dphiCOM;}
337 int thetaIDmin = m_thetaID[crysIDMax[0]];
338 if (m_thetaID[crysIDMax[1]] < m_thetaID[crysIDMax[0]]) {thetaIDmin = m_thetaID[crysIDMax[1]];}
339 getObjectPtr<TH2F>(
"dPhivsThetaID")->Fill(thetaIDmin + 0.001, dphiCOM);
340 if (dphiCOM<m_dPhiMin.at(thetaIDmin) or dphiCOM>m_dPhiMax.at(thetaIDmin)) {
return;}
345 memset(&EperCrys[0], 0, EperCrys.size()*
sizeof EperCrys[0]);
347 if (m_useCalDigits) {
348 for (
auto& eclCalDigit : m_eclCalDigitArray) {
349 int tempCrysID = eclCalDigit.getCellId() - 1;
350 EperCrys[tempCrysID] = eclCalDigit.getEnergy();
353 for (
auto& eclDigit : m_eclDigitArray) {
354 int tempCrysID = eclDigit.getCellId() - 1;
356 EperCrys[tempCrysID] = eclDigit.getAmp() * abs(ee5x5Calib[tempCrysID]) * ElectronicsCalib[tempCrysID];
362 for (
int ic = 0; ic < 2; ic++) {
363 int crysMax = crysIDMax[ic];
364 float expE = abs(Expee5x5E[crysMax]);
365 float sigmaExp = Expee5x5Sigma[crysMax];
366 std::vector<short int> neighbours = m_eclNeighbours5x5->getNeighbours(crysMax + 1);
369 double reducedExpE = expE;
371 for (
auto& cellID : neighbours) {
372 E25 += EperCrys[cellID - 1];
373 if (ee5x5Calib[cellID - 1] < 0.) {
374 reducedExpE -= EperCrys[cellID - 1];
379 double rexpE = reducedExpE / sigmaExp;
380 for (
auto& celli : neighbours) {
381 if (ee5x5Calib[celli - 1] > 0.) {
382 float rEi = EperCrys[celli - 1] / sigmaExp;
383 getObjectPtr<TH1F>(
"RvsCrysID")->Fill(celli - 0.999, rexpE * rEi);
384 getObjectPtr<TH1F>(
"NRvsCrysID")->Fill(celli - 0.999);
385 for (
auto& cellj : neighbours) {
386 if (ee5x5Calib[cellj - 1] > 0.) {
387 float rEj = EperCrys[cellj - 1] / sigmaExp;
388 getObjectPtr<TH2F>(
"Qmatrix")->Fill(celli - 0.999, cellj - 0.999, rEi * rEj);
395 getObjectPtr<TH2F>(
"EnVsCrysID")->Fill(crysMax + 0.001, E25 / expE);
396 getObjectPtr<TH1F>(
"EntriesvsCrys")->Fill(crysMax + 0.001);
397 getObjectPtr<TH1F>(
"ExpEvsCrys")->Fill(crysMax + 0.001, Expee5x5E[crysMax]);
398 getObjectPtr<TH1F>(
"ElecCalibvsCrys")->Fill(crysMax + 0.001, ElectronicsCalib[crysMax]);
399 getObjectPtr<TH1F>(
"InitialCalibvsCrys")->Fill(crysMax + 0.001, ee5x5Calib[crysMax]);
400 getObjectPtr<TH1F>(
"CalibEntriesvsCrys")->Fill(crysMax + 0.001);
Calibration collector module base class.
Class to provide momentum-related information from ECLClusters.
const TLorentzVector Get4MomentumFromCluster(const ECLCluster *cluster, ECLCluster::EHypothesisBit hypo)
Returns four momentum vector.
const TVector3 GetIPPosition()
Returns default IP position from beam parameters.
Class to store calibrated ECLDigits: ECLCalDigits.
int getCellId() const
Get Cell ID.
@ c_nPhotons
CR is split into n photons (N1)
Class to get the neighbours for a given cell id.
Calibration collector module that uses e+e- --> e+e- to do ECL single crystal energy calibration.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.