CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
DQA_TOF.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2#include "GaudiKernel/SmartDataPtr.h"
3#include "GaudiKernel/PropertyMgr.h"
4#include "GaudiKernel/Bootstrap.h"
5
6#include "GaudiKernel/INTupleSvc.h"
7#include "GaudiKernel/NTuple.h"
8#include "GaudiKernel/ITHistSvc.h"
9
10#include "CLHEP/Vector/ThreeVector.h"
11#include "CLHEP/Vector/LorentzVector.h"
12
14#include "EventModel/Event.h"
15
21#include "DQAEvent/DQAEvent.h"
22#include "DQA_TOF/DQA_TOF.h"
23#include <vector>
24/////////////////////////////////////////////////////////////////////////////
25DQA_TOF::DQA_TOF(const std::string& name, ISvcLocator* pSvcLocator) :
26 Algorithm(name, pSvcLocator) {
27
28 //Declare the properties
29 }
30// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
32 MsgStream log(msgSvc(), name());
33
34 log << MSG::INFO << "in initialize()" << endmsg;
35 StatusCode status;
36
37 StatusCode sc=service("THistSvc", m_thsvc);
38 if(sc.isFailure()) {
39 log << MSG::ERROR << "Couldn't get THistSvc" << endreq;
40 return StatusCode::FAILURE;
41 }
42
43 h_path=new TH1F("h_path","barrel ",200,0,200);
44 if(m_thsvc->regHist("/DQAHist/TOF/h_path",h_path).isFailure()){
45 log << MSG::ERROR << "Couldn't register h_path" << endreq;
46 }
47 h_Bzrhit=new TH1F("h_Bzrhit","barrel z hitmap",240,-120,120);
48 m_thsvc->regHist("/DQAHist/TOF/h_Bzrhit", h_Bzrhit);
49 h_Ezrhit=new TH1F("h_Ezrhit","endcap z hitmap",50,40,90);
50 m_thsvc->regHist("/DQAHist/TOF/h_Ezrhit", h_Ezrhit);
51 h_ph=new TH1F("h_ph","barrel Q",900,0,9000);
52 m_thsvc->regHist("/DQAHist/TOF/h_ph", h_ph);
53
54 W_delT=new TH2F("W_delT","barrel west PMT delT",176,0,176,300,-1.5,1.5);
55 m_thsvc->regHist("/DQAHist/TOF/W_delT", W_delT);
56 E_delT=new TH2F("E_delT","barrel east PMT delT",176,0,176,300,-1.5,1.5);
57 m_thsvc->regHist("/DQAHist/TOF/E_delT", E_delT);
58 counter=new TH2F("counter","barrel counter delT",176,0,176,300,-1.5,1.5);
59 m_thsvc->regHist("/DQAHist/TOF/counter", counter);
60 cluster=new TH2F("cluster","barrel cluster delT",88,0,88,300,-1.5,1.5);
61 m_thsvc->regHist("/DQAHist/TOF/cluster", cluster);
62 EC_delT=new TH2F("EC_delT","endcap delT",96,0,96,300,-1.5,1.5);
63 m_thsvc->regHist("/DQAHist/TOF/EC_delT", EC_delT);
64 Bt_delT=new TH1F("Bt_delT","barrel delT",300,-1.5,1.5);
65 m_thsvc->regHist("/DQAHist/TOF/Bt_delT", Bt_delT);
66 Et_delT=new TH1F("Et_delT","endcap delT",300,-1.5,1.5);
67 m_thsvc->regHist("/DQAHist/TOF/Et_delT", Et_delT);
68
69 B_path=new TH2F("B_path","barrel flight distance vs z",240,-120,120,200,0.0,200.0);
70 m_thsvc->regHist("/DQAHist/TOF/B_path", B_path);
71 E_path=new TH2F("E_path","endcap path distance vs z",50,40,90,200,0.0,200.0);
72 m_thsvc->regHist("/DQAHist/TOF/E_path", E_path);
73
74 delT_z1=new TH2F("delT_z1","barrel east delT vs Z",240,-120,120,300,-1.5,1.5);
75 m_thsvc->regHist("/DQAHist/TOF/delT_z1", delT_z1);
76 delT_z2=new TH2F("delT_z2","barrel west delT vs Z",240,-120,120,300,-1.5,1.5);
77 m_thsvc->regHist("/DQAHist/TOF/delT_z2", delT_z2);
78 delT_z3=new TH2F("delT_z3","barrel counter delT vs Z",240,-120,120,300,-1.5,1.5);
79 m_thsvc->regHist("/DQAHist/TOF/delT_z3", delT_z3);
80 delT_z4=new TH2F("delT_z4","barrel cluster delT vs Z",240,-120,120,300,-1.5,1.5);
81 m_thsvc->regHist("/DQAHist/TOF/delT_z4", delT_z4);
82
83 W_delT_Q=new TH2F("W_delT_Q","west barrel delT vs Q",900,0,9000,300,-1.5,1.5);
84 m_thsvc->regHist("/DQAHist/TOF/W_delT_Q", W_delT_Q);
85 E_delT_Q=new TH2F("E_delT_Q","east barrel delT vs Q",900,0,9000,300,-1.5,1.5);
86 m_thsvc->regHist("/DQAHist/TOF/E_delT_Q", E_delT_Q);
87
88 delT_pp=new TH1F("delT_pp","proton delT",300,-1.5,1.5);
89 m_thsvc->regHist("/DQAHist/TOF/delT_pp", delT_pp);
90 delT_pm=new TH1F("delT_Pm","anti-proton delT",300,-1.5,1.5);
91 m_thsvc->regHist("/DQAHist/TOF/delT_pm", delT_pm);
92 delT_pi=new TH1F("delT_pi","pi delT",300,-1.5,1.5);
93 m_thsvc->regHist("/DQAHist/TOF/delT_pi", delT_pi);
94 delT_k=new TH1F("delT_k","k delT",300,-1.5,1.5);
95 m_thsvc->regHist("/DQAHist/TOF/delT_k", delT_k);
96 log << MSG::INFO << "DQA_TOF successfully return from initialize()" <<endmsg;
97 return StatusCode::SUCCESS;
98}
99// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
100StatusCode DQA_TOF::execute() {
101 MsgStream log(msgSvc(), name());
102 log << MSG::INFO << "in execute()" << endreq;
103
104 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
105 int runNo=eventHeader->runNumber();
106 int event=eventHeader->eventNumber();
107 log<<MSG::DEBUG<<"run,evtnum="<<runNo<<","<<event<<endreq;
108
109 SmartDataPtr<DQAEvent::DQAEvent> dqaevt(eventSvc(), "/Event/DQATag");
110 if ( !dqaevt ) {
111 log << MSG::ERROR << "Error accessing DQAEvent" << endreq;
112 return StatusCode::FAILURE;
113 }
114
115 log << MSG::DEBUG << "event tag = " << dqaevt->EventTag() << endreq;
116
117 //if ( !dqaevt->Bhabha() )return StatusCode::SUCCESS;
118 SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc(), EventModel::EvtRec::EvtRecEvent);
119 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol(eventSvc(), EventModel::EvtRec::EvtRecTrackCol);
120 if(dqaevt->Bhabha()){
121 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
122 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
123 //if ( (*itTrk)->partId() != 1 ) continue; // only e+, e-
124 if(!(*itTrk)->isElectron())continue;
125 if(!(*itTrk)->isTofTrackValid())continue;
126 SmartRefVector<RecTofTrack> tofTrkCol = (*itTrk)->tofTrack();
127 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
128 for(;iter_tof != tofTrkCol.end(); iter_tof++ ) {
129 TofHitStatus *status = new TofHitStatus;
130 status->setStatus((*iter_tof)->status());
131 int tofid=(*iter_tof)->tofID();
132 double tof=(*iter_tof)->tof();
133 double texpE=(*iter_tof)->texpElectron();
134 double path=(*iter_tof)->path();
135 double Q=(*iter_tof)->ph();
136 double zrhit=(*iter_tof)->zrhit();
137
138 if(status->is_barrel()){
139 h_Bzrhit->Fill(zrhit);
140 h_path->Fill(path);
141 //h_ph->Fill(Q);
142 if(status->is_readout()){
143 h_ph->Fill(Q);
144 B_path->Fill(zrhit,path);
145 if(status->is_east()){ //barrel readout east
146 E_delT->Fill(tofid,tof-texpE);
147 delT_z1->Fill(zrhit,tof-texpE);
148 E_delT_Q->Fill(Q,tof-texpE);
149 }
150 else { //barrel readout west
151 W_delT->Fill(tofid,tof-texpE);
152 delT_z2->Fill(zrhit,tof-texpE);
153 W_delT_Q->Fill(Q,tof-texpE);
154 }
155 } //end readout
156 if(!status->is_readout()&&status->is_counter()){ //barrel counter
157 counter->Fill(tofid,tof-texpE);
158 delT_z3->Fill(zrhit,tof-texpE);
159 }
160 if(!status->is_counter()&&status->is_cluster()){ //barrel cluster
161 cluster->Fill(tofid,tof-texpE);
162 Bt_delT->Fill(tof-texpE);
163 delT_z4->Fill(zrhit,tof-texpE);
164 }
165 }
166 else{ //endcap
167 E_path->Fill(zrhit,path);
168 h_Ezrhit->Fill(zrhit);
169 EC_delT->Fill(tofid,tof-texpE);
170 Et_delT->Fill(tof-texpE);
171 }
172 } // loop 7 info fo tof
173 } // loop every track
174 }
175 else{
176 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
177 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
178 //if ( (*itTrk)->partId() != 1 ) continue; // only e+, e-
179 if(!(*itTrk)->isTofTrackValid())continue;
180 SmartRefVector<RecTofTrack> tofTrkCol = (*itTrk)->tofTrack();
181 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
182 for(;iter_tof != tofTrkCol.end(); iter_tof++ ) {
183 TofHitStatus *status = new TofHitStatus;
184 status->setStatus((*iter_tof)->status());
185 if(!status->is_cluster())continue;
186 double tof=(*iter_tof)->tof();
187 double texpPi=(*iter_tof)->texpPion();
188 double texpK=(*iter_tof)->texpKaon();
189 double texpP=(*iter_tof)->texpProton();
190 //if((*itTrk)->partId()==3){ //pi
191 if((*itTrk)->isPion()){
192 double texpPi=(*iter_tof)->texpPion();
193 delT_pi->Fill(tof-texpPi);
194 }
195 //else if((*itTrk)->partId()==4){ //k
196 else if((*itTrk)->isKaon()){
197 double texpK=(*iter_tof)->texpKaon();
198 delT_k->Fill(tof-texpK);
199 }
200 //else if((*itTrk)->partId()==5){ //p
201 else if((*itTrk)->isProton()){
202 double texpP=(*iter_tof)->texpProton();
203 if(!(*itTrk)->isMdcTrackValid())continue;
204 RecMdcTrack *mdcTrk = (*itTrk)->mdcTrack();
205 if(mdcTrk->charge()>0)delT_pp->Fill(tof-texpP);
206 else delT_pm->Fill(tof-texpP);
207 }
208 }
209 }
210 }
211 return StatusCode::SUCCESS;
212
213}
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
215StatusCode DQA_TOF::finalize() {
216 MsgStream log(msgSvc(), name());
217 log << MSG::INFO << "in finalize()" << endmsg;
218
219 return StatusCode::SUCCESS;
220}
221
222
int runNo
EvtRecTrackCol::iterator EvtRecTrackIterator
Definition: EvtRecTrack.h:131
IMessageSvc * msgSvc()
StatusCode finalize()
Definition: DQA_TOF.cxx:215
StatusCode execute()
Definition: DQA_TOF.cxx:100
DQA_TOF(const std::string &name, ISvcLocator *pSvcLocator)
Definition: DQA_TOF.cxx:25
StatusCode initialize()
Definition: DQA_TOF.cxx:31
const int charge() const
Definition: DstMdcTrack.h:53
bool is_barrel() const
Definition: TofHitStatus.h:26
bool is_cluster() const
Definition: TofHitStatus.h:25
void setStatus(unsigned int status)
bool is_counter() const
Definition: TofHitStatus.h:24
bool is_east() const
Definition: TofHitStatus.h:27
bool is_readout() const
Definition: TofHitStatus.h:23
_EXTERN_ std::string EvtRecEvent
Definition: EventModel.h:134
_EXTERN_ std::string EvtRecTrackCol
Definition: EventModel.h:135