BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibEvent Class Reference

#include <DedxCalibEvent.h>

+ Inheritance diagram for DedxCalibEvent:

Public Member Functions

 DedxCalibEvent (const std::string &name, ISvcLocator *pSvcLocator)
 
 ~DedxCalibEvent ()
 
void initializing ()
 
void BookHists ()
 
void genNtuple ()
 
void FillHists ()
 
void AnalyseHists ()
 
void WriteHists ()
 
- Public Member Functions inherited from DedxCalib
 DedxCalib (const std::string &name, ISvcLocator *pSvcLocator)
 
 ~DedxCalib ()
 
StatusCode initialize ()
 
StatusCode execute ()
 
StatusCode finalize ()
 
virtual void genNtuple ()=0
 
virtual void initializing ()=0
 
virtual void BookHists ()=0
 
virtual void FillHists ()=0
 
virtual void AnalyseHists ()=0
 
virtual void WriteHists ()=0
 

Public Attributes

float cut_wire
 
int m_count
 
int m_gap
 

Additional Inherited Members

- Protected Member Functions inherited from DedxCalib
double cal_dedx_bitrunc (float truncate, std::vector< double > phlist)
 
double cal_dedx (float truncate, std::vector< double > phlist)
 
void getCurvePar ()
 
void set_dEdx (int landau, float dEdx, int trkalg, int runflag, int dedxhit_use, float ptrk, float theta, float pl_rp, int vflag[3], double t0)
 
void ReadRecFileList ()
 
- Protected Attributes inherited from DedxCalib
IMdcGeomSvcgeosvc
 
IDedxCorrecSvcexsvc
 
float truncate
 
vector< double > Curve_Para
 
vector< double > Sigma_Para
 
int vFlag [3]
 
double m_dedx_exp [5]
 
double m_ex_sigma [5]
 
double m_pid_prob [5]
 
double m_chi_ex [5]
 
vector< string > m_recFileVector
 
int ParticleFlag
 
int m_calibflag
 
int m_phShapeFlag
 
std::string m_eventType
 
std::string m_recFileList
 
std::string m_rootfile
 
std::string m_curvefile
 

Detailed Description

Definition at line 8 of file DedxCalibEvent.h.

Constructor & Destructor Documentation

◆ DedxCalibEvent()

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

Definition at line 32 of file DedxCalibEvent.cxx.

32 : DedxCalib(name, pSvcLocator) {
33 //Declare the properties
34 declareProperty("CutWire", cut_wire=-1);
35 declareProperty("Count", m_count=1000000);
36 declareProperty("Gap", m_gap=1);
37}

◆ ~DedxCalibEvent()

DedxCalibEvent::~DedxCalibEvent ( )
inline

Definition at line 12 of file DedxCalibEvent.h.

12{};

Member Function Documentation

◆ AnalyseHists()

void DedxCalibEvent::AnalyseHists ( )
inlinevirtual

Implements DedxCalib.

Definition at line 20 of file DedxCalibEvent.h.

20{}

◆ BookHists()

void DedxCalibEvent::BookHists ( )
inlinevirtual

Implements DedxCalib.

Definition at line 17 of file DedxCalibEvent.h.

17{}

◆ FillHists()

void DedxCalibEvent::FillHists ( )
inlinevirtual

Implements DedxCalib.

Definition at line 19 of file DedxCalibEvent.h.

19{}

◆ genNtuple()

void DedxCalibEvent::genNtuple ( )
virtual

Implements DedxCalib.

Definition at line 131 of file DedxCalibEvent.cxx.

132{
133 MsgStream log(msgSvc(), name());
134 log << MSG::INFO << "DedxCalibEvent::genNtuple()" << endreq;
135
136 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
137 if (!eventHeader)
138 {
139 log << MSG::INFO << "Could not find Event Header" << endreq;
140 return;
141 }
142 int eventNO = eventHeader->eventNumber();
143 int runNO = eventHeader->runNumber();
144 // only save part of the events, in each cluster only m_count events saved, then jump with a gap
145 //std::cout << "eventNO " << eventNO << " evt_count " << evt_count << std::endl;
146 if(eventNO < evt_threshold) return;
147 evt_count ++;
148 if(evt_count == m_count){
149 evt_threshold = eventNO + m_gap;
150 evt_count = 0;
151 }
152 //std::cout << "evt_threshold " << evt_threshold << " evt_count " << evt_count << std::endl;
153
154
155 log << MSG::INFO << "now begin the event: runNO= "<<runNO<<" eventNO= "<< eventNO<< endreq;
156
157 int runFlag=0; //data type flag
158 if(runNO<RUN0) runFlag =0;
159 if(runNO>=RUN1 && runNO<RUN2) runFlag =1;
160 if(runNO>=RUN2 && runNO<RUN3) runFlag =2;
161 if(runNO>=RUN3 && runNO<RUN4) runFlag =3;
162 if(runNO>=RUN4 && runNO<RUN5) runFlag =4; //jpsi
163 if(runNO>=RUN5 && runNO<RUN6) runFlag =5; //psipp
164 if(runNO>=RUN6 && runNO<RUN7) runFlag =6; //psi4040, psip, jpsi...
165 if(runNO>=RUN7) runFlag =7; //psip
166
167 double tes = -999.0;
168 int esTimeflag = -1;
169 SmartDataPtr<RecEsTimeCol> aevtimeCol(eventSvc(),"/Event/Recon/RecEsTimeCol");
170 if( (!aevtimeCol) || (aevtimeCol->size()==0) ){
171 tes = -9999.0;
172 esTimeflag = -1;
173 }else{
174 RecEsTimeCol::iterator iter_evt = aevtimeCol->begin();
175 for(; iter_evt!=aevtimeCol->end(); iter_evt++){
176 tes = (*iter_evt)->getTest();
177 esTimeflag = (*iter_evt)->getStat();
178 }
179 }
180 if(runFlag ==2) {if( tes>1000 ) return;}
181 else if(runFlag ==3 ){if (tes>700 ) return;}
182 else {if (tes>1400 ) return;}
183
184 SmartDataPtr<EvtRecEvent>evtRecEvent(eventSvc(),"/Event/EvtRec/EvtRecEvent");
185 if(!evtRecEvent){
186 log << MSG::ERROR << "EvtRecEvent not found" << endreq;
187 return ;
188 }
189
190 SmartDataPtr<EvtRecTrackCol>
191 evtRecTrkCol(eventSvc(), "/Event/EvtRec/EvtRecTrackCol");
192 if(!evtRecTrkCol){
193 log << MSG::ERROR << "EvtRecTrackCol" << endreq;
194 return ;
195 }
196
197
198 Vint iGood;
199 iGood.clear();
200 int nCharge = 0;
201 double db=0,dz=0,pt0=0,charge0=0;
202 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
203 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
204 RecMdcDedx *it = (*itTrk)->mdcDedx();
205 if(it==0) continue;
206 HepVector a;
207 // if((*itTrk)->isMdcKalTrackValid()) cout << "mdckaltrack is valid" << endl;
208 // else cout << "mdckaltrack is not valid" << endl;
209
210 if((*itTrk)->isMdcKalTrackValid())
211 {
212 RecMdcKalTrack* trk = it->getMdcKalTrack();
213 if(ParticleFlag>-1&&ParticleFlag<5)
214 {
215 DstMdcKalTrack* dstTrk = trk;
216 a = dstTrk->getZHelix(ParticleFlag);
217 }
218 else
219 a = trk->getZHelix();
220 }
221 else if((*itTrk)->isMdcTrackValid())
222 {
223 RecMdcTrack* trk = it->getMdcTrack();
224 a = trk->helix();
225 }
226 else continue;
227 db = a[0];
228 dz = a[3];
229 pt0 = fabs(1.0/a[2]);
230 charge0 = ( a[2] > 0 )? 1 : -1;
231
232 //cout<<"db: "<<db<<" dz: "<<dz<<" pt0: "<<pt0<<" charge0: "<<charge0<<endl;
233 if(fabs(dz) >= VZ0CUT) continue;
234 if(db >= VR0CUT) continue;
235 if(pt0 >= PT0HighCut) continue;
236 if(pt0 <= PT0LowCut) continue;
237 iGood.push_back(it->trackId());
238 nCharge += charge0;
239 }
240
241
242 double poca_x=0,poca_y=0,poca_z=0;
243 float sintheta=0,costheta=0,ptrk=0,ptrk_t=0,charge=0,trackId=0;
244 int Nhit=0,usedhit=0,Nhits=0,Nphlisthit=0;
245 double dEdx_meas_hit=0, dEdx_meas=0,dEdx_meas_esat=0,dEdx_meas_norun=0;
246 double dEdx_hit[100]={0},pathlength_hit[100]={0},wid_hit[100]={0},layid_hit[100]={0},dd_in_hit[100]={0},eangle_hit[100]={0},zhit_hit[100]={0};
247 int trk_recalg = -1;
248 Identifier mdcid;
249
250 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
251 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
252 RecMdcDedx *it = (*itTrk)->mdcDedx();
253 if(it==0) continue;
254 bool flag = false;
255 for(unsigned int i = 0; i < iGood.size(); i++)
256 {
257 if(it->trackId()==iGood[i]) flag=true;
258 }
259 if(flag==false) continue;
260
261 HepVector a;
262 HepSymMatrix tkErrM;
263 if((*itTrk)->isMdcKalTrackValid())
264 {
265 poca_x = it->getMdcKalTrack()->x(); //get poca, default pid is pion; change pid using setPidType();
266 poca_y = it->getMdcKalTrack()->y();
267 poca_z = it->getMdcKalTrack()->z();
268
269 RecMdcKalTrack* trk = it->getMdcKalTrack();
270 if(ParticleFlag>-1&&ParticleFlag<5)
271 {
272 DstMdcKalTrack* dstTrk = trk;
273 a = dstTrk->getFHelix(ParticleFlag);
274 tkErrM = dstTrk->getFError(ParticleFlag);
275 }
276 else
277 {
278 a = trk->getFHelix();
279 tkErrM = trk->getFError();
280 }
281 }
282 else if((*itTrk)->isMdcTrackValid())
283 {
284 poca_x = it->getMdcTrack()->x();
285 poca_y = it->getMdcTrack()->y();
286 poca_z = it->getMdcTrack()->z();
287
288 RecMdcTrack* trk = it->getMdcTrack();
289 a = trk->helix();
290 tkErrM = trk->err();
291 }
292 else continue;
293
294 sintheta = sin(M_PI_2 - atan(a[4]));
295 costheta = cos(M_PI_2 - atan(a[4]));
296 ptrk_t = 1.0/fabs( a[2] );
297 ptrk = ptrk_t*sqrt(1 + a[4]*a[4]);
298 charge = ( a[2] > 0 )? 1 : -1;
299
300 Nhit = it->numTotalHits(); //total hits on track used as sample;
301 Nhits = (it->getVecDedxHits()).size(); //dedx hits on this track, they are put in phlist if layid>3
302 usedhit = it->numGoodHits(); //hits after truncting phlist and used in cal dE/dx value;
303 trk_recalg = it->status();
304 trackId = it->trackId();
305
306 if(m_eventType == "isBhabha")
307 {
308 if(runFlag ==3 &&(ptrk>1.88 || ptrk<1.80)) continue;
309 if(runFlag ==4 &&(ptrk>1.72 || ptrk<1.45)) continue;
310 if(runFlag ==5 &&(ptrk>2.00 || ptrk<1.70)) continue;
311 if(runFlag ==6 &&(ptrk>1.90 || ptrk<1.00)) continue;
312 if(runFlag ==7 &&(ptrk>1.90 || ptrk<0.90)) continue;
313 if(abs(costheta)>0.9) continue;
314
315 if(Nhit<20) continue;
316 if(usedhit<6) continue;
317 }
318
319
320 int layid=0,localwid=0,w0id=0,wid=0,lr=0;
321 double p_hit=0,adc_raw=0,tdc_raw=0,zhit=0,driftd=0,driftD=0,driftT=0,dd_in=0,dd_ex=0,eangle=0;
322 double ph=0,pathlength=0,rphi_path=0;
323 long k=0;
324
325 DedxHitRefVec gothits = it->getVecDedxHits();
326 DedxHitRefVec::iterator it_gothit = gothits.begin();
327 for(;it_gothit!=gothits.end(); it_gothit++)
328 {
329 if((*it_gothit)->isMdcHitValid())
330 {
331 RecMdcHit* itor = (*it_gothit)->getMdcHit();
332 mdcid = itor->getMdcId();
333 layid = MdcID::layer(mdcid);
334 localwid = MdcID::wire(mdcid);
335 w0id = geosvc->Layer(layid)->Wirst();
336 wid = w0id + localwid;
337 adc_raw = itor->getAdc();
338 tdc_raw = itor->getTdc();
339 zhit = itor->getZhit();
340
341 lr = itor->getFlagLR();
342 if(lr == 2) cout<<"lr = "<<lr<<endl;
343 if(lr == 0 || lr == 2) driftD = itor->getDriftDistLeft();
344 else driftD = itor->getDriftDistRight();
345 driftd = abs(driftD);
346 dd_in = itor->getDoca();
347 dd_ex = itor->getDoca();
348 if(lr == 0 || lr == 2 ) {dd_in = -abs(dd_in);dd_ex = -abs(dd_ex);}
349 if(lr == 1 ) {dd_in = abs(dd_in);dd_ex = abs(dd_ex);}
350 driftT = itor->getDriftT();
351 eangle = itor->getEntra();
352 eangle = eangle/M_PI;
353 }
354 else if((*it_gothit)->isMdcKalHelixSegValid())
355 {
356 RecMdcKalHelixSeg* itor = (*it_gothit)->getMdcKalHelixSeg();
357 HepVector a_hit_in = itor->getHelixIncl();
358 p_hit = 1.0/fabs(a_hit_in(3))*sqrt(1+a_hit_in(5)*a_hit_in(5));
359
360 mdcid = itor->getMdcId();
361 layid = MdcID::layer(mdcid);
362 localwid = MdcID::wire(mdcid);
363 w0id = geosvc->Layer(layid)->Wirst();
364 wid = w0id + localwid;
365 adc_raw = itor->getAdc();
366 tdc_raw = itor->getTdc();
367 zhit = itor->getZhit();
368
369 lr = itor->getFlagLR();
370 if(lr == 2) cout<<"lr = "<<lr<<endl;
371 driftD = itor->getDD();
372 driftd = abs(driftD);
373 driftT = itor->getDT();
374 dd_in = itor->getDocaIncl(); //getDocaIncl() include fit unused hit
375 dd_ex = itor->getDocaExcl(); //getDocaExcl() exclude fit unused hit
376 if(lr==-1 || lr == 2) {dd_in = dd_in; dd_ex = dd_ex;}
377 else if(lr ==1) {dd_in = -dd_in; dd_ex = -dd_ex;}
378 eangle = itor->getEntra();
379 eangle = eangle/M_PI;
380 }
381 else continue;
382
383 pathlength=(*it_gothit)->getPathLength();
384 rphi_path=pathlength*sintheta;
385 ph = (*it_gothit)->getDedx();
386 if(layid>3)
387 {
388 dEdx_hit[k]=adc_raw;
389 pathlength_hit[k]=pathlength;
390 wid_hit[k]=wid;
391 layid_hit[k]=layid;
392 dd_in_hit[k]=dd_in;
393 eangle_hit[k]=eangle;
394 zhit_hit[k]=zhit;
395
396 k++;
397 }
398
399 //cout<<"begin to Fill Ntuple n102!!!!!!!!!"<<endl;
400 m_charge1 = charge;
401 m_t01 = tes;
402 m_driftT = driftT;
403 m_tdc_raw = tdc_raw;
404 m_phraw = adc_raw;
405 m_exraw = ph;
406 m_localwid = localwid;
407 m_wire = wid;
408 m_runNO1 = runNO;
409 m_evtNO1 = eventNO;
410 m_runFlag1 = runFlag;
411 m_doca_in = dd_in;
412 m_doca_ex = dd_ex;
413 m_driftdist = driftD;
414 m_eangle = eangle;
415 m_zhit = zhit;
416 m_costheta1 = costheta;
417 m_pathL = pathlength;
418 m_layer = layid;
419 m_ptrk1 = ptrk;
420 m_ptrk_hit = p_hit;
421 m_trackId1 = trackId;
422 StatusCode status;
423 if(cut_wire>0){ if(cut_wire == m_wire) status = m_nt2->write(); }
424 else status = m_nt2->write();
425 if ( status.isFailure() )
426 log << MSG::ERROR << "Cannot fill Ntuple n102" << endreq;
427 }
428
429 Nphlisthit = k; //dedx hits on this track, exclude the first 3 layers
430 dEdx_meas = it->probPH();
431 dEdx_meas_esat = it->getDedxEsat();
432 dEdx_meas_norun = it->getDedxNoRun();
433 dEdx_meas_hit = it->getDedxHit();
434
435 //cout<<"begin to Fill Ntuple n103!!!!!!!"<<endl;
436 m_poca_x = poca_x;
437 m_poca_y = poca_y;
438 m_poca_z = poca_z;
439 m_ptrk_t=ptrk_t;
440 m_ptrk=ptrk;
441 m_sintheta=sintheta;
442 m_costheta=costheta;
443 m_charge=charge;
444 m_runNO = runNO;
445 m_runFlag = runFlag;
446 m_evtNO = eventNO;
447 m_t0 = tes;
448 m_trackId = trackId;
449 m_recalg = trk_recalg;
450
451 m_nhit=Nhit;
452 m_nhits=Nhits;
453 m_nphlisthit=Nphlisthit;
454 m_usedhit=usedhit;
455 for(int j=0; j<Nphlisthit; j++)
456 {
457 m_dEdx_hit[j]=dEdx_hit[j];
458 m_pathlength_hit[j]=pathlength_hit[j];
459 m_wid_hit[j]=wid_hit[j];
460 m_layid_hit[j]=layid_hit[j];
461 m_dd_in_hit[j]=dd_in_hit[j];
462 m_eangle_hit[j]=eangle_hit[j];
463 m_zhit_hit[j]=zhit_hit[j];
464 }
465
466 //m_dEdx_meas_hit = dEdx_meas_hit;
467 m_dEdx_meas = dEdx_meas;
468 //m_dEdx_meas_esat = dEdx_meas_esat;
469 //m_dEdx_meas_norun = dEdx_meas_norun;
470
471 m_parttype = it->particleId();
472 m_chidedxe=it->chiE();
473 m_chidedxmu=it->chiMu();
474 m_chidedxpi=it->chiPi();
475 m_chidedxk=it->chiK();
476 m_chidedxp=it->chiP();
477 for(int i=0;i<5;i++)
478 {
479 m_probpid[i]= it->getPidProb(i);
480 m_expectid[i]= it->getDedxExpect(i);
481 m_sigmaid[i]= it->getSigmaDedx(i);
482 }
483 StatusCode status;
484 if(cut_wire<0) status = m_nt1->write();
485 if ( status.isFailure() )
486 {
487 log << MSG::ERROR << "Cannot fill Ntuple n103" << endreq;
488 }
489 }
490 //cout<<"track iteration ended!!!!!!!!!!"<<endl;
491}
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
static int evt_threshold(0)
const long int RUN5
const double VZ0CUT
const double PT0HighCut
const long int RUN6
const long int RUN2
const long int RUN7
const long int RUN1
const long int RUN4
const double VR0CUT
const double PT0LowCut
const long int RUN3
const long int RUN0
EvtRecTrackCol::iterator EvtRecTrackIterator
Definition: EvtRecTrack.h:111
std::vector< int > Vint
Definition: Gam4pikp.cxx:52
SmartRefVector< RecMdcDedxHit > DedxHitRefVec
Definition: RecMdcDedx.h:27
IMessageSvc * msgSvc()
#define M_PI
Definition: TConstant.h:4
int ParticleFlag
Definition: DedxCalib.h:50
IMdcGeomSvc * geosvc
Definition: DedxCalib.h:29
std::string m_eventType
Definition: DedxCalib.h:53
double probPH() const
Definition: DstMdcDedx.h:66
double chiE() const
Definition: DstMdcDedx.h:59
int particleId() const
Definition: DstMdcDedx.h:33
int numTotalHits() const
Definition: DstMdcDedx.h:65
int numGoodHits() const
Definition: DstMdcDedx.h:64
int status() const
Definition: DstMdcDedx.h:56
double chiPi() const
Definition: DstMdcDedx.h:61
double chiK() const
Definition: DstMdcDedx.h:62
double chiMu() const
Definition: DstMdcDedx.h:60
int trackId() const
Definition: DstMdcDedx.h:32
double chiP() const
Definition: DstMdcDedx.h:63
const double y() const
const HepVector & getZHelix(const int pid) const
const HepSymMatrix & getFError(const int pid) const
const HepVector & getFHelix(const int pid) const
const double z() const
const double x() const
const HepSymMatrix err() const
const HepVector helix() const
......
const double z() const
Definition: DstMdcTrack.h:63
const double y() const
Definition: DstMdcTrack.h:62
const double x() const
Definition: DstMdcTrack.h:61
virtual const MdcGeoLayer *const Layer(unsigned id)=0
int Wirst(void) const
Definition: MdcGeoLayer.h:157
static int layer(const Identifier &id)
Values of different levels (failure returns 0)
Definition: MdcID.cxx:49
static int wire(const Identifier &id)
Definition: MdcID.cxx:54
RecMdcKalTrack * getMdcKalTrack(void)
Definition: RecMdcDedx.h:73
double getDedxNoRun(void)
Definition: RecMdcDedx.h:65
DedxHitRefVec getVecDedxHits() const
Definition: RecMdcDedx.h:61
double getPidProb(int pid) const
Definition: RecMdcDedx.h:70
RecMdcTrack * getMdcTrack(void)
Definition: RecMdcDedx.h:72
double getDedxEsat(void)
Definition: RecMdcDedx.h:64
double getSigmaDedx(int pid) const
Definition: RecMdcDedx.h:69
double getDedxExpect(int pid) const
Definition: RecMdcDedx.h:68
double getDedxHit(void)
Definition: RecMdcDedx.h:63
const int getFlagLR(void) const
Definition: RecMdcHit.h:47
const double getZhit(void) const
Definition: RecMdcHit.h:55
const double getAdc(void) const
Definition: RecMdcHit.h:51
const Identifier getMdcId(void) const
Definition: RecMdcHit.h:49
const double getTdc(void) const
Definition: RecMdcHit.h:50
const double getDriftDistRight(void) const
Definition: RecMdcHit.h:43
const double getDriftT(void) const
Definition: RecMdcHit.h:52
const double getEntra(void) const
Definition: RecMdcHit.h:54
const double getDriftDistLeft(void) const
Definition: RecMdcHit.h:42
const double getDoca(void) const
Definition: RecMdcHit.h:53
Identifier getMdcId(void) const
int getFlagLR(void) const
HepVector & getHelixIncl(void)
double getEntra(void) const
double getDD(void) const
double getDocaIncl(void) const
double getDocaExcl(void) const
double getTdc(void) const
double getZhit(void) const
double getDT(void) const
double getAdc(void) const
const HepVector & getZHelix() const
const HepSymMatrix & getFError() const
const HepVector & getFHelix() const
float costheta
float charge
float dEdx_meas
float ptrk
float dEdx_hit[100]

◆ initializing()

void DedxCalibEvent::initializing ( )
virtual

Implements DedxCalib.

Definition at line 40 of file DedxCalibEvent.cxx.

41{
42 MsgStream log(msgSvc(), name());
43 log << MSG::INFO << "DedxCalibEvent::initializing()" << endreq;
44
45 StatusCode status;
46 NTuplePtr nt1(ntupleSvc(),"FILE100/n103");
47 if ( nt1 )
48 m_nt1 = nt1;
49 else
50 {
51 m_nt1= ntupleSvc()->book("FILE100/n103",CLID_ColumnWiseTuple,"dEdx per track");
52 if ( m_nt1 )
53 {
54 status = m_nt1->addItem("ptrk",m_ptrk);
55 status = m_nt1->addItem("ptrk_t",m_ptrk_t);
56 status = m_nt1->addItem("sintheta",m_sintheta);
57 status = m_nt1->addItem("costheta",m_costheta);
58 status = m_nt1->addItem("charge",m_charge);
59 status = m_nt1->addItem("runNO",m_runNO);
60 status = m_nt1->addItem("runFlag",m_runFlag);
61 status = m_nt1->addItem("evtNO",m_evtNO);
62 status = m_nt1->addItem("t0",m_t0);
63 status = m_nt1->addItem("trackId",m_trackId);
64 status = m_nt1->addItem("poca_x",m_poca_x);
65 status = m_nt1->addItem("poca_y",m_poca_y);
66 status = m_nt1->addItem("poca_z",m_poca_z);
67 status = m_nt1->addItem("recalg",m_recalg);
68 status = m_nt1->addItem("nhit",m_nhit);
69 status = m_nt1->addItem("nhits",m_nhits);
70 status = m_nt1->addItem("usedhit",m_usedhit);
71
72 status = m_nt1->addItem("ndedxhit",m_nphlisthit,0,100);
73 status = m_nt1->addIndexedItem("dEdx_hit",m_nphlisthit,m_dEdx_hit);
74 status = m_nt1->addIndexedItem("pathlength_hit",m_nphlisthit,m_pathlength_hit);
75 status = m_nt1->addIndexedItem("wid_hit",m_nphlisthit,m_wid_hit);
76 status = m_nt1->addIndexedItem("layid_hit",m_nphlisthit,m_layid_hit);
77 status = m_nt1->addIndexedItem("dd_in_hit",m_nphlisthit,m_dd_in_hit);
78 status = m_nt1->addIndexedItem("eangle_hit",m_nphlisthit,m_eangle_hit);
79 status = m_nt1->addIndexedItem("zhit_hit",m_nphlisthit,m_zhit_hit);
80
81 //status = m_nt1->addItem("dEdx_meas_hit", m_dEdx_meas_hit);
82 status = m_nt1->addItem("dEdx_meas", m_dEdx_meas);
83 //status = m_nt1->addItem("dEdx_meas_esat", m_dEdx_meas_esat);
84 //status = m_nt1->addItem("dEdx_meas_norun", m_dEdx_meas_norun);
85
86 status = m_nt1->addItem("type",m_parttype);
87 status = m_nt1->addItem("chidedx_e",m_chidedxe);
88 status = m_nt1->addItem("chidedx_mu",m_chidedxmu);
89 status = m_nt1->addItem("chidedx_pi",m_chidedxpi);
90 status = m_nt1->addItem("chidedx_k",m_chidedxk);
91 status = m_nt1->addItem("chidedx_p",m_chidedxp);
92 status = m_nt1->addItem("partid",5,m_probpid);
93 status = m_nt1->addItem("expectid",5,m_expectid);
94 status = m_nt1->addItem("sigmaid",5,m_sigmaid);
95 }
96 }
97
98 NTuplePtr nt2(ntupleSvc(),"FILE100/n102");
99 if ( nt2 ) m_nt2 = nt2;
100 else
101 {
102 m_nt2= ntupleSvc()->book("FILE100/n102",CLID_RowWiseTuple,"dE/dx per hit");
103 if ( m_nt2 )
104 {
105 status = m_nt2->addItem("charge",m_charge1);
106 status = m_nt2->addItem("adc_raw",m_phraw);
107 status = m_nt2->addItem("exraw",m_exraw);
108 status = m_nt2->addItem("runNO",m_runNO1);
109 status = m_nt2->addItem("evtNO",m_evtNO1);
110 status = m_nt2->addItem("runFlag",m_runFlag1);
111 status = m_nt2->addItem("wire",m_wire);
112 status = m_nt2->addItem("doca_in",m_doca_in);
113 status = m_nt2->addItem("doca_ex",m_doca_ex);
114 status = m_nt2->addItem("driftdist",m_driftdist);
115 status = m_nt2->addItem("eangle",m_eangle);
116 status = m_nt2->addItem("zhit",m_zhit);
117 status = m_nt2->addItem("costheta1",m_costheta1);
118 status = m_nt2->addItem("path_rphi",m_pathL);
119 status = m_nt2->addItem("layer",m_layer);
120 status = m_nt2->addItem("ptrk1",m_ptrk1);
121 status = m_nt2->addItem("ptrk_hit",m_ptrk_hit);
122 status = m_nt2->addItem("t01",m_t01);
123 status = m_nt2->addItem("tdc_raw",m_tdc_raw);
124 status = m_nt2->addItem("driftT",m_driftT);
125 status = m_nt2->addItem("localwid",m_localwid);
126 status = m_nt2->addItem("trackId1",m_trackId1);
127 }
128 }
129}
INTupleSvc * ntupleSvc()

◆ WriteHists()

void DedxCalibEvent::WriteHists ( )
inlinevirtual

Implements DedxCalib.

Definition at line 21 of file DedxCalibEvent.h.

21{}

Member Data Documentation

◆ cut_wire

float DedxCalibEvent::cut_wire

Definition at line 13 of file DedxCalibEvent.h.

Referenced by DedxCalibEvent(), and genNtuple().

◆ m_count

int DedxCalibEvent::m_count

Definition at line 14 of file DedxCalibEvent.h.

Referenced by DedxCalibEvent(), and genNtuple().

◆ m_gap

int DedxCalibEvent::m_gap

Definition at line 15 of file DedxCalibEvent.h.

Referenced by DedxCalibEvent(), and genNtuple().


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