Belle II Software development
eclee5x5CollectorModule.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
9/* Own header. */
10#include <ecl/modules/eclee5x5Collector/eclee5x5CollectorModule.h>
11
12/* ECL headers. */
13#include <ecl/dataobjects/ECLCalDigit.h>
14#include <ecl/dataobjects/ECLDigit.h>
15#include <ecl/dataobjects/ECLElementNumbers.h>
16#include <ecl/dbobjects/ECLCrystalCalib.h>
17
18/* Basf2 headers. */
19#include <analysis/ClusterUtility/ClusterUtils.h>
20#include <framework/dataobjects/EventMetaData.h>
21#include <framework/datastore/RelationVector.h>
22#include <framework/geometry/VectorUtil.h>
23#include <mdst/dataobjects/ECLCluster.h>
24#include <mdst/dataobjects/TRGSummary.h>
25
26/* ROOT headers. */
27#include <Math/Vector3D.h>
28#include <Math/Vector4D.h>
29#include <TH2F.h>
30
31using namespace std;
32using namespace Belle2;
33
34
35//-----------------------------------------------------------------
36// Register the Module
37//-----------------------------------------------------------------
38REG_MODULE(eclee5x5Collector);
39
40//-----------------------------------------------------------------
41// Implementation
42//-----------------------------------------------------------------
43
44//-----------------------------------------------------------------------------------------------------
45
47 m_ECLExpee5x5E("ECLExpee5x5E"), m_ElectronicsCalib("ECLCrystalElectronics"), m_ee5x5Calib("ECLCrystalEnergyee5x5"),
48 m_selectdPhiData("ECLeedPhiData"), m_selectdPhiMC("ECLeedPhiMC")
49{
50 // Set module properties
51 setDescription("Calibration Collector Module for ECL single crystal energy calibration using Bhabha events");
53 addParam("thetaLabMinDeg", m_thetaLabMinDeg, "minimum 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);
62 addParam("expectedEnergyScale", m_expectedEnergyScale, "scale expected energies for non-4S calibration", 1.);
63}
64
65
66
70{
71
74 B2INFO("eclee5x5Collector: Experiment = " << m_evtMetaData->getExperiment() << " run = " << m_evtMetaData->getRun() << " sqrts = "
75 << m_sqrts);
76
77
80 auto EnVsCrysID = new TH2F("EnVsCrysID", "Normalized 5x5 energy for each crystal;crystal ID;E25/Expected",
82 1.2);
83 registerObject<TH2F>("EnVsCrysID", EnVsCrysID);
84
85 auto RvsCrysID = new TH1F("RvsCrysID", "E_exp x E_crysID / sigma^2;crysID;sum of E_exp x E_crysID/sigma^2",
87 registerObject("RvsCrysID", RvsCrysID);
88
89 auto NRvsCrysID = new TH1F("NRvsCrysID", "Entries in RvsCrysID vs crysID;crysID;Entries in RvsCrysID",
91 registerObject("NRvsCrysID", NRvsCrysID);
92
93 auto Qmatrix = new TH2F("Qmatrix", "E_i x E_j/sigma^2;crysID i;crysID j", ECLElementNumbers::c_NCrystals, 0,
95 registerObject("Qmatrix", Qmatrix);
96
97
98 auto ElecCalibvsCrys = new TH1F("ElecCalibvsCrys", "Sum electronics calib const vs crystal ID;crystal ID;calibration constant",
100 registerObject<TH1F>("ElecCalibvsCrys", ElecCalibvsCrys);
101 auto ExpEvsCrys = new TH1F("ExpEvsCrys", "Sum expected energy calib const vs crystalID;crystal ID;calibration constant",
104 registerObject<TH1F>("ExpEvsCrys", ExpEvsCrys);
105 auto InitialCalibvsCrys = new TH1F("InitialCalibvsCrys", "Sum initial calib const vs crystal ID;crystal ID;calibration constant",
107 registerObject<TH1F>("InitialCalibvsCrys", InitialCalibvsCrys);
108
109 auto CalibEntriesvsCrys = new TH1F("CalibEntriesvsCrys", "Entries in calib vs crys histograms;crystal ID;Entries per crystal",
112 registerObject<TH1F>("CalibEntriesvsCrys", CalibEntriesvsCrys);
113
114 auto EntriesvsCrys = new TH1F("EntriesvsCrys", "Selected Bhabha clusters vs crystal ID;crystal ID;Entries",
116 registerObject<TH1F>("EntriesvsCrys", EntriesvsCrys);
117
118 auto dPhivsThetaID = new TH2F("dPhivsThetaID",
119 "Phi* vs thetaID forward, pass thetaSum,E0,E1;thetaID of forward cluster;dPhi COM (deg)", 69, 0, 69, 150, 165, 180);
120 registerObject<TH2F>("dPhivsThetaID", dPhivsThetaID);
121
122 //------------------------------------------------------------------------
124 B2INFO("Input parameters to eclee5x5Collector:");
125 B2INFO("thetaLabMinDeg: " << m_thetaLabMinDeg);
126 B2INFO("thetaLabMaxDeg: " << m_thetaLabMaxDeg);
127 m_thetaLabMin = m_thetaLabMinDeg / TMath::RadToDeg();
128 m_thetaLabMax = m_thetaLabMaxDeg / TMath::RadToDeg();
129 B2INFO("minE0: " << m_minE0);
130 B2INFO("minE1: " << m_minE1);
131 B2INFO("maxdThetaSum: " << m_maxdThetaSum);
132 B2INFO("dPhiScale: " << m_dPhiScale);
133 B2INFO("maxTime: " << m_maxTime);
134 B2INFO("useCalDigits: " << m_useCalDigits);
135 B2INFO("requireL1: " << m_requireL1);
136 B2INFO("expectedEnergyScale: " << m_expectedEnergyScale);
137
141
144
147 if (m_ECLExpee5x5E.hasChanged()) {
148 Expee5x5E = m_ECLExpee5x5E->getCalibVector();
149 Expee5x5Sigma = m_ECLExpee5x5E->getCalibUncVector();
150 }
151 if (m_ElectronicsCalib.hasChanged()) {ElectronicsCalib = m_ElectronicsCalib->getCalibVector();}
152 if (m_ee5x5Calib.hasChanged()) {ee5x5Calib = m_ee5x5Calib->getCalibVector();}
153 if (m_selectdPhiMC.hasChanged() and m_useCalDigits) {
154 meandPhi = m_selectdPhiMC->getCalibVector();
155 widthdPhi = m_selectdPhiMC->getCalibUncVector();
156 } else if (m_selectdPhiData.hasChanged()) {
157 meandPhi = m_selectdPhiData->getCalibVector();
158 widthdPhi = m_selectdPhiData->getCalibUncVector();
159 }
160
162 for (int ic = 1; ic < 9000; ic += 1000) {
163 B2INFO("DB constants for cellID=" << ic << ": ee5x5Calib = " << ee5x5Calib[ic - 1] << " Expee5x5E = " << Expee5x5E[ic - 1] <<
164 " ElectronicsCalib = " <<
165 ElectronicsCalib[ic - 1]);
166 }
167
169 for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
170 if (ElectronicsCalib[crysID] <= 0) {B2FATAL("eclee5x5Collector: ElectronicsCalib = " << ElectronicsCalib[crysID] << " for crysID = " << crysID);}
171 if (Expee5x5E[crysID] == 0) {B2FATAL("eclee5x5Collector: Expee5x5E = 0 for crysID = " << crysID);}
172 if (ee5x5Calib[crysID] == 0) {B2FATAL("eclee5x5Collector: ee5x5Calib = 0 for crysID = " << crysID);}
173 }
174
177 m_eclClusterArray.isRequired();
178 m_eclCalDigitArray.isRequired();
179 if (!m_useCalDigits) {m_eclDigitArray.isRequired();}
180 m_evtMetaData.isRequired();
181
183 //..Derive ThetaID of each crystal, and cut on dPhi* as a function of thetaID
184 int crysID = 0;
185 for (int it = 0; it < 69; it++) {
186
187 //..dPhi* cuts are actually a function of thetaID, not crysID
188 m_dPhiMin.push_back(meandPhi[crysID] - m_dPhiScale * widthdPhi[crysID]);
189 m_dPhiMax.push_back(meandPhi[crysID] + m_dPhiScale * widthdPhi[crysID]);
190 for (int ic = 0; ic < m_eclNeighbours5x5->getCrystalsPerRing(it); ic++) {
191 m_thetaID.at(crysID) = it;
192 crysID++;
193 }
194 }
195
196}
197
198
202{
203
205 if (storeCalib) {
206 for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
207 getObjectPtr<TH1F>("ExpEvsCrys")->Fill(crysID + 0.001, Expee5x5E[crysID]);
208 getObjectPtr<TH1F>("ElecCalibvsCrys")->Fill(crysID + 0.001, ElectronicsCalib[crysID]);
209 getObjectPtr<TH1F>("InitialCalibvsCrys")->Fill(crysID + 0.001, ee5x5Calib[crysID]);
210 getObjectPtr<TH1F>("CalibEntriesvsCrys")->Fill(crysID + 0.001);
211 }
212 storeCalib = false;
213 }
214
217 bool newConst = false;
218 if (m_ECLExpee5x5E.hasChanged()) {
219 newConst = true;
220 B2INFO("ECLExpee5x5E has changed, exp = " << m_evtMetaData->getExperiment() << " run = " << m_evtMetaData->getRun());
221 Expee5x5E = m_ECLExpee5x5E->getCalibVector();
222 Expee5x5Sigma = m_ECLExpee5x5E->getCalibUncVector();
223 }
224 if (m_ElectronicsCalib.hasChanged()) {
225 newConst = true;
226 B2INFO("ECLCrystalElectronics has changed, exp = " << m_evtMetaData->getExperiment() << " run = " << m_evtMetaData->getRun());
227 ElectronicsCalib = m_ElectronicsCalib->getCalibVector();
228 }
229 if (m_ee5x5Calib.hasChanged()) {
230 newConst = true;
231 B2INFO("ECLCrystalEnergyee5x5 has changed, exp = " << m_evtMetaData->getExperiment() << " run = " << m_evtMetaData->getRun());
232 ee5x5Calib = m_ee5x5Calib->getCalibVector();
233 }
234
235 if (newConst) {
236 for (int ic = 1; ic < 9000; ic += 1000) {
237 B2INFO("DB constants for cellID=" << ic << ": ee5x5Calib = " << ee5x5Calib[ic - 1] << " Expee5x5E = " << Expee5x5E[ic - 1] <<
238 " ElectronicsCalib = " <<
239 ElectronicsCalib[ic - 1]);
240 }
241
243 for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
244 if (ElectronicsCalib[crysID] <= 0) {B2FATAL("eclee5x5Collector: ElectronicsCalib = " << ElectronicsCalib[crysID] << " for crysID = " << crysID);}
245 if (Expee5x5E[crysID] == 0) {B2FATAL("eclee5x5Collector: Expee5x5E = 0 for crysID = " << crysID);}
246 if (ee5x5Calib[crysID] == 0) {B2FATAL("eclee5x5Collector: ee5x5Calib = 0 for crysID = " << crysID);}
247 }
248 }
249
252 if (m_requireL1) {
253 unsigned int L1TriggerResults = m_TRGResults->getTRGSummary(0);
254 if (L1TriggerResults == 0) {return;}
255 }
256
257 //------------------------------------------------------------------------
259 int icMax[2] = { -1, -1};
260 double maxClustE[2] = { -1., -1.};
261 int nclust = m_eclClusterArray.getEntries();
262 for (int ic = 0; ic < nclust; ic++) {
265 if (eClust > maxClustE[0]) {
266 maxClustE[1] = maxClustE[0];
267 icMax[1] = icMax[0];
268 maxClustE[0] = eClust;
269 icMax[0] = ic;
270 } else if (eClust > maxClustE[1]) {
271 maxClustE[1] = eClust;
272 icMax[1] = ic;
273 }
274 }
275 }
276
277 //------------------------------------------------------------------------
280 if (icMax[1] == -1) {return;}
281 double theta0 = m_eclClusterArray[icMax[0]]->getTheta();
282 double theta1 = m_eclClusterArray[icMax[1]]->getTheta();
283 if (theta0 < m_thetaLabMin || theta0 > m_thetaLabMax || theta1 < m_thetaLabMin || theta1 > m_thetaLabMax) {return;}
284
286 double t0 = m_eclClusterArray[icMax[0]]->getTime();
287 double dt990 = m_eclClusterArray[icMax[0]]->getDeltaTime99();
288 double t1 = m_eclClusterArray[icMax[1]]->getTime();
289 double dt991 = m_eclClusterArray[icMax[1]]->getDeltaTime99();
290 double dt99min = dt990;
291 if (dt991 < dt990) {dt99min = dt991;}
292 if (dt99min <= 0) {dt99min = 0.0001;}
293 if (abs(t1 - t0) > dt99min * m_maxTime) {return;}
294
295 //------------------------------------------------------------------------
297 ClusterUtils cUtil;
298 const ROOT::Math::XYZVector clustervertex = cUtil.GetIPPosition();
299
300 double phi0 = m_eclClusterArray[icMax[0]]->getPhi();
301 ROOT::Math::XYZVector p30;
302 VectorUtil::setMagThetaPhi(p30, maxClustE[0], theta0, phi0);
303 const ROOT::Math::PxPyPzEVector p40 = cUtil.Get4MomentumFromCluster(m_eclClusterArray[icMax[0]], clustervertex,
305
306 double phi1 = m_eclClusterArray[icMax[1]]->getPhi();
307 ROOT::Math::XYZVector p31;
308 VectorUtil::setMagThetaPhi(p31, maxClustE[1], theta1, phi1);
309 const ROOT::Math::PxPyPzEVector p41 = cUtil.Get4MomentumFromCluster(m_eclClusterArray[icMax[1]], clustervertex,
311
313 ROOT::Math::PxPyPzEVector p40COM = m_boostrotate.rotateLabToCms() * p40;
314 ROOT::Math::PxPyPzEVector p41COM = m_boostrotate.rotateLabToCms() * p41;
315 double theta01COM = (p41COM.Theta() + p40COM.Theta()) * TMath::RadToDeg();
316 if (abs(theta01COM - 180.) > m_maxdThetaSum) {return;}
317
319 if (p40COM.E() < m_minE0 * m_sqrts and p41COM.E() < m_minE0 * m_sqrts) {return;}
320 if (p40COM.E() < m_minE1 * m_sqrts or p41COM.E() < m_minE1 * m_sqrts) {return;}
321
322
323 //------------------------------------------------------------------------
325 int crysIDMax[2] = { -1, -1};
326 double crysEMax[2] = { -1., -1.};
327 for (int imax = 0; imax < 2; imax++) {
328 auto eclClusterRelations = m_eclClusterArray[icMax[imax]]->getRelationsTo<ECLCalDigit>("ECLCalDigits");
329 for (unsigned int ir = 0; ir < eclClusterRelations.size(); ir++) {
330 const auto calDigit = eclClusterRelations.object(ir);
331 int tempCrysID = calDigit->getCellId() - 1;
332 float tempE = calDigit->getEnergy();
333 if (tempE > crysEMax[imax]) {
334 crysEMax[imax] = tempE;
335 crysIDMax[imax] = tempCrysID;
336 }
337 }
338 }
339
341 double dphiCOM = abs(p41COM.Phi() - p40COM.Phi()) * TMath::RadToDeg();
342 if (dphiCOM > 180.) {dphiCOM = 360. - dphiCOM;}
343
344 int thetaIDmin = m_thetaID[crysIDMax[0]];
345 if (m_thetaID[crysIDMax[1]] < m_thetaID[crysIDMax[0]]) {thetaIDmin = m_thetaID[crysIDMax[1]];}
346 getObjectPtr<TH2F>("dPhivsThetaID")->Fill(thetaIDmin + 0.001, dphiCOM);
347 if (dphiCOM<m_dPhiMin.at(thetaIDmin) or dphiCOM>m_dPhiMax.at(thetaIDmin)) {return;}
348
349
350 //------------------------------------------------------------------------
352 memset(&EperCrys[0], 0, EperCrys.size()*sizeof EperCrys[0]);
353
354 if (m_useCalDigits) {
355 for (auto& eclCalDigit : m_eclCalDigitArray) {
356 int tempCrysID = eclCalDigit.getCellId() - 1;
357 EperCrys[tempCrysID] = eclCalDigit.getEnergy();
358 }
359 } else {
360 for (auto& eclDigit : m_eclDigitArray) {
361 int tempCrysID = eclDigit.getCellId() - 1;
363 EperCrys[tempCrysID] = eclDigit.getAmp() * abs(ee5x5Calib[tempCrysID]) * ElectronicsCalib[tempCrysID];
364 }
365 }
366
367 //------------------------------------------------------------------------
368 //** Quantities needed for the 5x5 calibration */
369 for (int ic = 0; ic < 2; ic++) {
370 int crysMax = crysIDMax[ic];
371 float expE = m_expectedEnergyScale * abs(Expee5x5E[crysMax]);
372 float sigmaExp = m_expectedEnergyScale * Expee5x5Sigma[crysMax];
373 std::vector<short int> neighbours = m_eclNeighbours5x5->getNeighbours(crysMax + 1);
374
375 //** Energy in 5x5, and expected energy corrected for crystals that will not be calibrated */
376 double reducedExpE = expE;
377 double E25 = 0.;
378 for (auto& cellID : neighbours) {
379 E25 += EperCrys[cellID - 1];
380 if (ee5x5Calib[cellID - 1] < 0.) {
381 reducedExpE -= EperCrys[cellID - 1];
382 }
383 }
384
385 //** now the vector and matrix used in the calibration */
386 double rexpE = reducedExpE / sigmaExp;
387 for (auto& celli : neighbours) {
388 if (ee5x5Calib[celli - 1] > 0.) {
389 float rEi = EperCrys[celli - 1] / sigmaExp;
390 getObjectPtr<TH1F>("RvsCrysID")->Fill(celli - 0.999, rexpE * rEi);
391 getObjectPtr<TH1F>("NRvsCrysID")->Fill(celli - 0.999);
392 for (auto& cellj : neighbours) {
393 if (ee5x5Calib[cellj - 1] > 0.) {
394 float rEj = EperCrys[cellj - 1] / sigmaExp;
395 getObjectPtr<TH2F>("Qmatrix")->Fill(celli - 0.999, cellj - 0.999, rEi * rEj);
396 }
397 }
398 }
399 }
400
401 //** Record normalized energy and corresponding calib values. Expee5x5E is negative if the algorithm was unable to calculate a value. In this case, the nominal input value has been stored with a minus sign */
402 getObjectPtr<TH2F>("EnVsCrysID")->Fill(crysMax + 0.001, E25 / expE);
403 getObjectPtr<TH1F>("EntriesvsCrys")->Fill(crysMax + 0.001);
404 getObjectPtr<TH1F>("ExpEvsCrys")->Fill(crysMax + 0.001, Expee5x5E[crysMax]);
405 getObjectPtr<TH1F>("ElecCalibvsCrys")->Fill(crysMax + 0.001, ElectronicsCalib[crysMax]);
406 getObjectPtr<TH1F>("InitialCalibvsCrys")->Fill(crysMax + 0.001, ee5x5Calib[crysMax]);
407 getObjectPtr<TH1F>("CalibEntriesvsCrys")->Fill(crysMax + 0.001);
408 }
409}
Calibration collector module base class.
void registerObject(std::string name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
Class to provide momentum-related information from ECLClusters.
Definition: ClusterUtils.h:36
const ROOT::Math::PxPyPzEVector Get4MomentumFromCluster(const ECLCluster *cluster, ECLCluster::EHypothesisBit hypo)
Returns four momentum vector.
Definition: ClusterUtils.cc:25
const ROOT::Math::XYZVector GetIPPosition()
Returns default IP position from beam parameters.
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
int getCellId() const
Get Cell ID.
Definition: ECLCalDigit.h:118
@ c_nPhotons
CR is split into n photons (N1)
Class to get the neighbours for a given cell id.
Definition: ECLNeighbours.h:25
const std::vector< short int > & getNeighbours(short int cid) const
Return the neighbours for a given cell ID.
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Definition: ECLNeighbours.h:39
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
double getCMSEnergy() const
Returns CMS energy of e+e- (aka.
const ROOT::Math::LorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
DBObjPtr< ECLCrystalCalib > m_selectdPhiData
dPhi cut
StoreArray< ECLDigit > m_eclDigitArray
Required input array of ECLDigits.
std::vector< float > widthdPhi
width of requirement on dPhi from DB object
std::vector< float > Expee5x5E
vector of energies obtained from DB object
bool m_requireL1
require events to satisfy a level 1 trigger (false)
DBObjPtr< ECLCrystalCalib > m_ECLExpee5x5E
Expected energies from database.
double m_maxTime
maximum cluster time diff abs(t1-t0)/dt99 (10)
eclee5x5CollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
double m_maxdThetaSum
abs(theta0* + theta1* - 180 deg) must be less than less (2 deg)
StoreArray< ECLCluster > m_eclClusterArray
Required arrays.
StoreObjPtr< TRGSummary > m_TRGResults
dataStore TRGSummary
bool storeCalib
force the input calibration constants to be saved first event
double m_minE0
minimum energy of cluster 0: E*0/sqrts (0.45)
std::vector< float > EperCrys
Energy for each crystal from ECLDigit or ECLCalDigit (GeV)
std::vector< int > m_thetaID
thetaID of each crystal
ECL::ECLNeighbours * m_eclNeighbours5x5
Neighbour map of 25 crystals.
bool m_useCalDigits
use eclCalDigit to determine MC deposited energy (false)
std::vector< float > Expee5x5Sigma
vector of sigmaE obtained from DB object
std::vector< float > m_dPhiMin
minimum dPhi* as a function of thetaID
std::vector< float > m_dPhiMax
maximum dPhi* as a function of thetaID
void collect() override
Select events and crystals and accumulate histograms.
StoreObjPtr< EventMetaData > m_evtMetaData
dataStore EventMetaData
std::vector< float > ElectronicsCalib
vector obtained from DB object
StoreArray< ECLCalDigit > m_eclCalDigitArray
Required input array of ECLCalDigits.
double m_thetaLabMax
m_thetaLabMaxDeg converted to radians
std::vector< float > meandPhi
mean requirement on dPhi from DB object
double m_sqrts
sqrt s from m_boostrotate
void prepare() override
Define histograms and read payloads from DB.
std::vector< float > ee5x5Calib
vector obtained from DB object
double m_thetaLabMin
Some other useful quantities.
DBObjPtr< ECLCrystalCalib > m_selectdPhiMC
DB object for MC.
double m_thetaLabMaxDeg
maximum ecl cluster theta in lab (150 degrees)
DBObjPtr< ECLCrystalCalib > m_ee5x5Calib
Existing single crystal calibration from DB; will be updated by CAF.
double m_thetaLabMinDeg
Parameters to control the job.
PCmsLabTransform m_boostrotate
boost from COM to lab and visa versa
double m_dPhiScale
scale dPhi* cut by this factor (1)
double m_minE1
minimum energy of cluster 1: E*1/sqrts (0.40)
DBObjPtr< ECLCrystalCalib > m_ElectronicsCalib
Electronics calibration from database.
double m_expectedEnergyScale
scale expected energies for non-4S calibration (1.)
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
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
STL namespace.