Belle II Software prerelease-10-00-00a
CDCDedxDQM.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#include <cdc/modules/CDCDedxDQM/CDCDedxDQM.h>
10
11#include <TDirectory.h>
12
13using namespace Belle2;
14
15REG_MODULE(CDCDedxDQM);
16
17//---------------------------------
19{
20 setPropertyFlags(c_ParallelProcessingCertified); // parallel processing
21 setDescription("CDC dE/dx DQM plots with bhabha/hadron events.");
22 addParam("mmode", mmode, "default monitoring mode is basic", std::string("basic"));
23}
24
25//---------------------------------
27{
28
29 TDirectory* oldDir = gDirectory;
30 oldDir->mkdir("CDCDedx");
31 oldDir->cd("CDCDedx");
32
33 int expNum = -1;
34 int runNum = -1;
35 double rungain = -99.0;
36
37 if (m_MetaDataPtr) {
38 expNum = int(m_MetaDataPtr->getExperiment());
39 runNum = int(m_MetaDataPtr->getRun());
40 if (m_DBRunGain) rungain = m_DBRunGain->getRunGain();
41 }
42
43 hMeta = new TH1D("hMeta", "hMeta", 3, 0.5, 3.5);
44 hMeta->GetXaxis()->SetTitle("Quantity");
45 hMeta->GetYaxis()->SetTitle("Values");
46 hMeta->SetTitle(Form("(Exp:%d, Run:%d, RG:%0.03f)", expNum, runNum, rungain));
47 hMeta->GetXaxis()->SetBinLabel(1, "nevt");
48 hMeta->GetXaxis()->SetBinLabel(2, "nbhabha");
49 hMeta->GetXaxis()->SetBinLabel(3, "nhadron");
50
51 hdEdx = new TH1D("hdEdx", ";CDC dE/dx;Entries", 100, 0., 2.5);
52 hinjtimeHer = new TH2D("hinjtimeHer", ";injection time (#mu s); CDC dE/dx", 40, 0, 80e3, 50, 0, 2.5);
53 hinjtimeLer = new TH2D("hinjtimeLer", ";injection time (#mu s); CDC dE/dx", 40, 0, 80e3, 50, 0, 2.5);
54 hdEdxvsP = new TH2D("hdEdxVsP", ";#it{p}_{CDC} (GeV/c);CDC dE/dx", 100, 0.05, 2.50, 100, 0.35, 10.0);
55 hdEdxvsEvt = new TH2D("hdEdxvsEvt", ";Events(M);CDC dE/dx", 50, 0, 200, 50, 0.00, 2.0);
56 hdEdxvsCosth = new TH2D("hdEdxvsCosth", ";cos#theta (e^{-}e^{+} tracks);CDC dE/dx", 50, -1.00, 1.00, 50, 0.00, 2.5);
57 hdEdxvsPhi = new TH2D("hdEdxvsPhi", ";#phi (e^{-}e^{+} tracks);CDC dE/dx", 50, -3.20, 3.20, 50, 0.00, 2.5);
58 if (mmode != "basic") {
59 hWires = new TH2F("hWires", "All Wires;", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
60 hWires->GetXaxis()->SetTitle("CDC-wire map: counter-clockwise and start from +x");
61 hWireStatus = new TH2F("hWireStatus", "Wire Status", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
62 hWireStatus->GetXaxis()->SetTitle("CDC-wire map: counter-clockwise and start from +x");
63 }
64 oldDir->cd();
65
66}
67
68
69//---------------------------------
71{
72
73 if (!m_cdcDedxTracks.isOptional()) {
74 B2WARNING("Missing CDCDedxTracks array, CDCDedxDQM is skipped.");
75 return;
76 }
77
78 m_TrgResult.isOptional();
79 m_cdcDedxTracks.isRequired();
80 REG_HISTOGRAM
81
82}
83
84//-------------------------------
86{
87
88 if (!m_cdcDedxTracks.isOptional()) {
89 B2WARNING("Missing CDCDedxTracks array, CDCDedxDQM is skipped.");
90 return;
91 }
92
93 hMeta->Reset();
94 hdEdx->Reset();
95 hinjtimeHer->Reset();
96 hinjtimeLer->Reset();
97 hdEdxvsP->Reset();
98 hdEdxvsCosth->Reset();
99 hdEdxvsPhi->Reset();
100 hdEdxvsEvt->Reset();
101 if (mmode != "basic") {
102 hWires->Reset();
103 hWireStatus->Reset();
104 }
105}
106
107
108//----------------------------
110{
111
112 if (!m_cdcDedxTracks.isOptional()) return;
113
114 if (!m_TrgResult.isValid()) {
115 B2WARNING("Required SoftwareTriggerResult object not available: CDCDedxDQM is skipped");
116 return;
117 }
118
119 const std::map<std::string, int>& fresults = m_TrgResult->getResults();
120 if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()
121 and fresults.find("software_trigger_cut&skim&accept_hadron") == fresults.end())return;
122
123 const bool IsBhabhaEvt = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
125 const bool IsHadronEvt = (m_TrgResult->getResult("software_trigger_cut&skim&accept_hadron") ==
127
128 m_nEvt += 1;
129 if (!IsBhabhaEvt and !IsHadronEvt)return;
130 if (IsBhabhaEvt)m_nBEvt += 1;
131 if (IsHadronEvt)m_nHEvt += 1;
132
133 //Get current evt number
134 int event = -1;
135
136 if (m_MetaDataPtr)event = int(m_MetaDataPtr->getEvent());
137
138 for (int idedx = 0; idedx < m_cdcDedxTracks.getEntries(); idedx++) {
139
140 CDCDedxTrack* dedxTrack = m_cdcDedxTracks[idedx];
141 if (!dedxTrack || dedxTrack->size() == 0)continue;
142
143 const Track* track = dedxTrack->getRelatedFrom<Track>();
144 if (!track)continue;
145
146 const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
147 if (!fitResult)continue;
148
149 UncertainHelix helix = fitResult->getUncertainHelix();
150 static DBObjPtr<BeamSpot> beamSpotDB;
151 helix.passiveMoveBy(ROOT::Math::XYZVector(beamSpotDB->getIPPosition()));
152 const auto& frame = ReferenceFrame::GetCurrent();
153 double dr = frame.getVertex(ROOT::Math::XYZVector(helix.getPerigee())).Rho();
154 double dz = frame.getVertex(ROOT::Math::XYZVector(helix.getPerigee())).Z();
155 if (dr >= 1.0 || fabs(dz) >= 1.0)continue;
156
157 //CDC acceptance as well
158 double costh = dedxTrack->getCosTheta();
159 if (costh < TMath::Cos(150.0 * TMath::DegToRad()))continue;
160 if (costh > TMath::Cos(17.0 * TMath::DegToRad())) continue;
161
162 //clean tracks with relaxed nhit cut
163 double nhits = dedxTrack->getNLayerHits();
164 if (costh > -0.55 && costh < 0.820) {
165 if (nhits < 20)continue;
166 } else {
167 if (costh <= -0.62 || costh >= 0.880) {
168 if (nhits < 8)continue;
169 if (costh > 0 && nhits < 10)continue;
170 } else {
171 if (nhits < 15)continue;
172 }
173 }
174
175 double dedxnosat = dedxTrack->getDedxNoSat();
176 if (dedxnosat < 0)continue;
177
178 double dedx = dedxTrack->getDedx();
179 if (dedx < 0)continue;
180
181 double pCDC = dedxTrack->getMomentum();
182 if (pCDC <= 0) continue;
183
184 double pTrk = fitResult->getMomentum().R();
185 if (pTrk <= 0) continue;
186
187 if (IsBhabhaEvt) {
188 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
189 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
190 double TrkEoverP = eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons) / pTrk;
191 if (TrkEoverP > 0) {
192 if (abs(TrkEoverP - 1.0) > 0.25)continue;
193 }
194 }
195
196 hdEdx->Fill(dedxnosat);
197 double phi = fitResult->getMomentum().Phi();
198 if (hdEdxvsCosth->Integral() <= 80000)hdEdxvsCosth->Fill(costh, dedxnosat);
199 if (hdEdxvsPhi->Integral() <= 80000)hdEdxvsPhi->Fill(phi, dedxnosat);
200
201 if (event >= 150e6)event = 150e6 - 100;
202 event = int(event / 5e5);
203 hdEdxvsEvt->Fill(event, dedxnosat);
204
205 // And check if the stored data is valid and if an injection happened recently
206 if (TTDInfo->isValid() && TTDInfo->hasInjection()) {
207 if (TTDInfo->isHER())
208 hinjtimeHer->Fill(TTDInfo->getTimeSinceLastInjectionInMicroSeconds(), dedxnosat);
209 else
210 hinjtimeLer->Fill(TTDInfo->getTimeSinceLastInjectionInMicroSeconds(), dedxnosat);
211 } else
212 return;
213
214 }
215 if (IsHadronEvt && hdEdxvsP->Integral() <= 80000)hdEdxvsP->Fill(pCDC, dedx);
216
217 if (mmode != "basic") {
218 for (int ihit = 0; ihit < dedxTrack->size(); ++ihit) {
219 int iwire = dedxTrack->getWire(ihit);
220 double iadc = dedxTrack->getADCCount(ihit);
221 if (m_adc[iwire].size() < 50)m_adc[iwire].push_back(iadc); //just contiung dead
222 }
223 }
224 }
225
226}
227
228//---------------------------------
230{
231
232 hMeta->SetBinContent(1, m_nEvt);
233 hMeta->SetBinContent(2, m_nBEvt);
234 hMeta->SetBinContent(3, m_nHEvt);
235
236 if (hdEdx->GetEntries() > 0) {
237 hdEdx->GetXaxis()->SetRange(hdEdx->FindFirstBinAbove(0, 1), hdEdx->FindLastBinAbove(0, 1));
238 }
239
240 if (hdEdxvsEvt->GetEntries() > 0) {
241 hdEdxvsEvt->GetXaxis()->SetRange(hdEdxvsEvt->FindFirstBinAbove(0, 1), hdEdxvsEvt->FindLastBinAbove(0, 1));
242 }
243
244 if (hinjtimeHer->GetEntries() > 0) {
245 hinjtimeHer->GetXaxis()->SetRange(hinjtimeHer->FindFirstBinAbove(0, 1), hinjtimeHer->FindLastBinAbove(0, 1));
246 }
247 if (hinjtimeLer->GetEntries() > 0) {
248 hinjtimeLer->GetXaxis()->SetRange(hinjtimeLer->FindFirstBinAbove(0, 1), hinjtimeLer->FindLastBinAbove(0, 1));
249 }
250 //get dead wire pattern
251 if (mmode != "basic") plotWireMap();
252}
253
254
255//---------------------------------
257{
258 B2INFO("CDCDedxDQMModule: terminate called");
259}
260
261//------------------------------------
263{
264
265 B2INFO("Creating CDCGeometryPar object");
267
268 int jwire = -1;
269 int nbadwires = 0;
270
271 for (unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
272 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
273 jwire++;
274 double phi = 2.*TMath::Pi() * (iwire / double(cdcgeo.nWiresInLayer(ilay)));
275 double radius = cdcgeo.senseWireR(ilay) / 100.;
276 double x = radius * cos(phi);
277 double y = radius * sin(phi);
278 hWires->Fill(x, y);
279 if (m_adc[jwire].size() > 0)continue;
280 nbadwires++;
281 hWireStatus->Fill(x, y);
282 }
283 }
284 hWireStatus->SetTitle(Form("%d", nbadwires));
285}
TH2F * hWireStatus
dead wire status
Definition CDCDedxDQM.h:106
DBObjPtr< CDCDedxRunGain > m_DBRunGain
Run gain DB object.
Definition CDCDedxDQM.h:109
TH2D * hdEdxvsCosth
dedx vs costh
Definition CDCDedxDQM.h:104
virtual void initialize() override
Initialize the module.
Definition CDCDedxDQM.cc:70
TH2D * hdEdxvsP
dedx vs p
Definition CDCDedxDQM.h:102
TH2F * hWires
all wire mapping
Definition CDCDedxDQM.h:107
virtual void event() override
This method is called for each event.
virtual void endRun() override
This method is called at the end of each run.
virtual void terminate() override
End of the event processing.
TH1D * hMeta
metadata
Definition CDCDedxDQM.h:98
StoreObjPtr< EventMetaData > m_MetaDataPtr
Store array for metadata info.
Definition CDCDedxDQM.h:85
virtual void beginRun() override
This method is called for each run.
Definition CDCDedxDQM.cc:85
StoreArray< CDCDedxTrack > m_cdcDedxTracks
Store array for CDCDedxTrack.
Definition CDCDedxDQM.h:87
int m_nEvt
accepted events
Definition CDCDedxDQM.h:90
TH2D * hdEdxvsPhi
dedx vs phi
Definition CDCDedxDQM.h:103
int m_nBEvt
bhabha events
Definition CDCDedxDQM.h:91
std::array< std::vector< double >, c_nSenseWires > m_adc
adc per wire for wire status
Definition CDCDedxDQM.h:94
TH2D * hdEdxvsEvt
dedx vs event
Definition CDCDedxDQM.h:105
int m_nHEvt
hadron events
Definition CDCDedxDQM.h:92
StoreObjPtr< EventLevelTriggerTimeInfo > TTDInfo
Store array for injection time info.
Definition CDCDedxDQM.h:88
CDCDedxDQMModule()
Default constructor.
Definition CDCDedxDQM.cc:18
TH2D * hinjtimeLer
injection time in LER
Definition CDCDedxDQM.h:101
TH2D * hinjtimeHer
injection time in HER
Definition CDCDedxDQM.h:100
void plotWireMap()
function to plot wire status map (all, bad)
std::string mmode
monitoring mode all/basic
Definition CDCDedxDQM.h:96
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
Definition CDCDedxDQM.h:86
virtual void defineHisto() override
Definition of histograms.
Definition CDCDedxDQM.cc:26
Debug output for CDCDedxPID module.
int getADCCount(int i) const
Return the adcCount for this hit.
double getDedx() const
Get dE/dx truncated mean for this track.
int getNLayerHits() const
Return the number of layer hits for this track.
double getCosTheta() const
Return cos(theta) for this track.
int getWire(int i) const
Return the sensor ID for this hit: wire number for CDC (0-14336)
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
int size() const
Return the number of hits for this track.
double getMomentum() const
Return the track momentum valid in the CDC.
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
static const ChargedStable pion
charged pion particle
Definition Const.h:661
Class for accessing objects in the database.
Definition DBObjPtr.h:21
ECL cluster data.
Definition ECLCluster.h:27
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition ECLCluster.h:351
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition ECLCluster.cc:23
@ c_nPhotons
CR is split into n photons (N1)
Definition ECLCluster.h:41
HistoModule()
Constructor.
Definition HistoModule.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
static const ReferenceFrame & GetCurrent()
Get current rest frame.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Values of the result of a track fit with a given particle hypothesis.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
UncertainHelix getUncertainHelix() const
Conversion to framework Uncertain Helix (i.e., with covariance).
Class that bundles various TrackFitResults.
Definition Track.h:25
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
double passiveMoveBy(const ROOT::Math::XYZVector &by)
Moves origin of the coordinate system (passive transformation) by the given vector.
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.