BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
DQA_TOF Class Reference

#include <DQA_TOF.h>

+ Inheritance diagram for DQA_TOF:

Public Member Functions

 DQA_TOF (const std::string &name, ISvcLocator *pSvcLocator)
 
StatusCode initialize ()
 
StatusCode execute ()
 
StatusCode finalize ()
 

Detailed Description

Definition at line 8 of file DQA_TOF.h.

Constructor & Destructor Documentation

◆ DQA_TOF()

DQA_TOF::DQA_TOF ( const std::string &  name,
ISvcLocator *  pSvcLocator 
)

Definition at line 26 of file DQA_TOF.cxx.

26 :
27 Algorithm(name, pSvcLocator) {
28
29 //Declare the properties
30 }

Member Function Documentation

◆ execute()

StatusCode DQA_TOF::execute ( )

Definition at line 101 of file DQA_TOF.cxx.

101 {
102 MsgStream log(msgSvc(), name());
103 log << MSG::INFO << "in execute()" << endreq;
104
105 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
106 int runNo=eventHeader->runNumber();
107 int event=eventHeader->eventNumber();
108 log<<MSG::DEBUG<<"run,evtnum="<<runNo<<","<<event<<endreq;
109
110 SmartDataPtr<DQAEvent::DQAEvent> dqaevt(eventSvc(), "/Event/DQATag");
111 if ( !dqaevt ) {
112 log << MSG::ERROR << "Error accessing DQAEvent" << endreq;
113 return StatusCode::FAILURE;
114 }
115
116 log << MSG::DEBUG << "event tag = " << dqaevt->EventTag() << endreq;
117
118 //if ( !dqaevt->Bhabha() )return StatusCode::SUCCESS;
119 SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc(), EventModel::EvtRec::EvtRecEvent);
120 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol(eventSvc(), EventModel::EvtRec::EvtRecTrackCol);
121 if(dqaevt->Bhabha()){
122 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
123 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
124 //if ( (*itTrk)->partId() != 1 ) continue; // only e+, e-
125 if(!(*itTrk)->isElectron())continue;
126 if(!(*itTrk)->isTofTrackValid())continue;
127 SmartRefVector<RecTofTrack> tofTrkCol = (*itTrk)->tofTrack();
128 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
129 for(;iter_tof != tofTrkCol.end(); iter_tof++ ) {
130 TofHitStatus *status = new TofHitStatus;
131 status->setStatus((*iter_tof)->status());
132 int tofid=(*iter_tof)->tofID();
133 double tof=(*iter_tof)->tof();
134 double texpE=(*iter_tof)->texpElectron();
135 double path=(*iter_tof)->path();
136 double Q=(*iter_tof)->ph();
137 double zrhit=(*iter_tof)->zrhit();
138
139 if(status->is_barrel()){
140 h_Bzrhit->Fill(zrhit);
141 h_path->Fill(path);
142 //h_ph->Fill(Q);
143 if(status->is_readout()){
144 h_ph->Fill(Q);
145 B_path->Fill(zrhit,path);
146 if(status->is_east()){ //barrel readout east
147 E_delT->Fill(tofid,tof-texpE);
148 delT_z1->Fill(zrhit,tof-texpE);
149 E_delT_Q->Fill(Q,tof-texpE);
150 }
151 else { //barrel readout west
152 W_delT->Fill(tofid,tof-texpE);
153 delT_z2->Fill(zrhit,tof-texpE);
154 W_delT_Q->Fill(Q,tof-texpE);
155 }
156 } //end readout
157 if(!status->is_readout()&&status->is_counter()){ //barrel counter
158 counter->Fill(tofid,tof-texpE);
159 delT_z3->Fill(zrhit,tof-texpE);
160 }
161 if(!status->is_counter()&&status->is_cluster()){ //barrel cluster
162 cluster->Fill(tofid,tof-texpE);
163 Bt_delT->Fill(tof-texpE);
164 delT_z4->Fill(zrhit,tof-texpE);
165 }
166 }
167 else{ //endcap
168 E_path->Fill(zrhit,path);
169 h_Ezrhit->Fill(zrhit);
170 EC_delT->Fill(tofid,tof-texpE);
171 Et_delT->Fill(tof-texpE);
172 }
173 } // loop 7 info fo tof
174 } // loop every track
175 }
176 else{
177 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
178 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
179 //if ( (*itTrk)->partId() != 1 ) continue; // only e+, e-
180 if(!(*itTrk)->isTofTrackValid())continue;
181 SmartRefVector<RecTofTrack> tofTrkCol = (*itTrk)->tofTrack();
182 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
183 for(;iter_tof != tofTrkCol.end(); iter_tof++ ) {
184 TofHitStatus *status = new TofHitStatus;
185 status->setStatus((*iter_tof)->status());
186 if(!status->is_cluster())continue;
187 double tof=(*iter_tof)->tof();
188 double texpPi=(*iter_tof)->texpPion();
189 double texpK=(*iter_tof)->texpKaon();
190 double texpP=(*iter_tof)->texpProton();
191 //if((*itTrk)->partId()==3){ //pi
192 if((*itTrk)->isPion()){
193 double texpPi=(*iter_tof)->texpPion();
194 delT_pi->Fill(tof-texpPi);
195 }
196 //else if((*itTrk)->partId()==4){ //k
197 else if((*itTrk)->isKaon()){
198 double texpK=(*iter_tof)->texpKaon();
199 delT_k->Fill(tof-texpK);
200 }
201 //else if((*itTrk)->partId()==5){ //p
202 else if((*itTrk)->isProton()){
203 double texpP=(*iter_tof)->texpProton();
204 if(!(*itTrk)->isMdcTrackValid())continue;
205 RecMdcTrack *mdcTrk = (*itTrk)->mdcTrack();
206 if(mdcTrk->charge()>0)delT_pp->Fill(tof-texpP);
207 else delT_pm->Fill(tof-texpP);
208 }
209 }
210 }
211 }
212 return StatusCode::SUCCESS;
213
214}
int runNo
Definition: DQA_TO_DB.cxx:12
EvtRecTrackCol::iterator EvtRecTrackIterator
Definition: EvtRecTrack.h:111
IMessageSvc * msgSvc()
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:116
_EXTERN_ std::string EvtRecTrackCol
Definition: EventModel.h:117

◆ finalize()

StatusCode DQA_TOF::finalize ( )

Definition at line 216 of file DQA_TOF.cxx.

216 {
217 MsgStream log(msgSvc(), name());
218 log << MSG::INFO << "in finalize()" << endmsg;
219
220 return StatusCode::SUCCESS;
221}

◆ initialize()

StatusCode DQA_TOF::initialize ( )

Definition at line 32 of file DQA_TOF.cxx.

32 {
33 MsgStream log(msgSvc(), name());
34
35 log << MSG::INFO << "in initialize()" << endmsg;
36 StatusCode status;
37
38 StatusCode sc=service("THistSvc", m_thsvc);
39 if(sc.isFailure()) {
40 log << MSG::ERROR << "Couldn't get THistSvc" << endreq;
41 return StatusCode::FAILURE;
42 }
43
44 h_path=new TH1F("h_path","barrel ",200,0,200);
45 if(m_thsvc->regHist("/DQAHist/TOF/h_path",h_path).isFailure()){
46 log << MSG::ERROR << "Couldn't register h_path" << endreq;
47 }
48 h_Bzrhit=new TH1F("h_Bzrhit","barrel z hitmap",240,-120,120);
49 m_thsvc->regHist("/DQAHist/TOF/h_Bzrhit", h_Bzrhit);
50 h_Ezrhit=new TH1F("h_Ezrhit","endcap z hitmap",50,40,90);
51 m_thsvc->regHist("/DQAHist/TOF/h_Ezrhit", h_Ezrhit);
52 h_ph=new TH1F("h_ph","barrel Q",900,0,9000);
53 m_thsvc->regHist("/DQAHist/TOF/h_ph", h_ph);
54
55 W_delT=new TH2F("W_delT","barrel west PMT delT",176,0,176,300,-1.5,1.5);
56 m_thsvc->regHist("/DQAHist/TOF/W_delT", W_delT);
57 E_delT=new TH2F("E_delT","barrel east PMT delT",176,0,176,300,-1.5,1.5);
58 m_thsvc->regHist("/DQAHist/TOF/E_delT", E_delT);
59 counter=new TH2F("counter","barrel counter delT",176,0,176,300,-1.5,1.5);
60 m_thsvc->regHist("/DQAHist/TOF/counter", counter);
61 cluster=new TH2F("cluster","barrel cluster delT",88,0,88,300,-1.5,1.5);
62 m_thsvc->regHist("/DQAHist/TOF/cluster", cluster);
63 EC_delT=new TH2F("EC_delT","endcap delT",96,0,96,300,-1.5,1.5);
64 m_thsvc->regHist("/DQAHist/TOF/EC_delT", EC_delT);
65 Bt_delT=new TH1F("Bt_delT","barrel delT",300,-1.5,1.5);
66 m_thsvc->regHist("/DQAHist/TOF/Bt_delT", Bt_delT);
67 Et_delT=new TH1F("Et_delT","endcap delT",300,-1.5,1.5);
68 m_thsvc->regHist("/DQAHist/TOF/Et_delT", Et_delT);
69
70 B_path=new TH2F("B_path","barrel flight distance vs z",240,-120,120,200,0.0,200.0);
71 m_thsvc->regHist("/DQAHist/TOF/B_path", B_path);
72 E_path=new TH2F("E_path","endcap path distance vs z",50,40,90,200,0.0,200.0);
73 m_thsvc->regHist("/DQAHist/TOF/E_path", E_path);
74
75 delT_z1=new TH2F("delT_z1","barrel east delT vs Z",240,-120,120,300,-1.5,1.5);
76 m_thsvc->regHist("/DQAHist/TOF/delT_z1", delT_z1);
77 delT_z2=new TH2F("delT_z2","barrel west delT vs Z",240,-120,120,300,-1.5,1.5);
78 m_thsvc->regHist("/DQAHist/TOF/delT_z2", delT_z2);
79 delT_z3=new TH2F("delT_z3","barrel counter delT vs Z",240,-120,120,300,-1.5,1.5);
80 m_thsvc->regHist("/DQAHist/TOF/delT_z3", delT_z3);
81 delT_z4=new TH2F("delT_z4","barrel cluster delT vs Z",240,-120,120,300,-1.5,1.5);
82 m_thsvc->regHist("/DQAHist/TOF/delT_z4", delT_z4);
83
84 W_delT_Q=new TH2F("W_delT_Q","west barrel delT vs Q",900,0,9000,300,-1.5,1.5);
85 m_thsvc->regHist("/DQAHist/TOF/W_delT_Q", W_delT_Q);
86 E_delT_Q=new TH2F("E_delT_Q","east barrel delT vs Q",900,0,9000,300,-1.5,1.5);
87 m_thsvc->regHist("/DQAHist/TOF/E_delT_Q", E_delT_Q);
88
89 delT_pp=new TH1F("delT_pp","proton delT",300,-1.5,1.5);
90 m_thsvc->regHist("/DQAHist/TOF/delT_pp", delT_pp);
91 delT_pm=new TH1F("delT_Pm","anti-proton delT",300,-1.5,1.5);
92 m_thsvc->regHist("/DQAHist/TOF/delT_pm", delT_pm);
93 delT_pi=new TH1F("delT_pi","pi delT",300,-1.5,1.5);
94 m_thsvc->regHist("/DQAHist/TOF/delT_pi", delT_pi);
95 delT_k=new TH1F("delT_k","k delT",300,-1.5,1.5);
96 m_thsvc->regHist("/DQAHist/TOF/delT_k", delT_k);
97 log << MSG::INFO << "DQA_TOF successfully return from initialize()" <<endmsg;
98 return StatusCode::SUCCESS;
99}

The documentation for this class was generated from the following files: