BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
inclks.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2#include "GaudiKernel/AlgFactory.h"
3#include "GaudiKernel/ISvcLocator.h"
4#include "GaudiKernel/SmartDataPtr.h"
5#include "GaudiKernel/IDataProviderSvc.h"
6#include "GaudiKernel/PropertyMgr.h"
7#include "GaudiKernel/Bootstrap.h"
8
9#include "GaudiKernel/INTupleSvc.h"
10#include "GaudiKernel/NTuple.h"
11#include "GaudiKernel/ITHistSvc.h"
12//#include "GaudiKernel/IHistogramSvc.h"
13
15#include "EventModel/Event.h"
16
21
22#include "TMath.h"
23#include "CLHEP/Vector/ThreeVector.h"
24#include "CLHEP/Vector/LorentzVector.h"
25#include "CLHEP/Vector/TwoVector.h"
26
27using CLHEP::Hep3Vector;
28using CLHEP::Hep2Vector;
29using CLHEP::HepLorentzVector;
30#include "CLHEP/Geometry/Point3D.h"
31#ifndef ENABLE_BACKWARDS_COMPATIBILITY
33#endif
34
37#include "VertexFit/VertexFit.h"
40
41//
42#include "DQAinclKs/inclks.h"
43
44#include <vector>
45
46const double mk0 = 0.497648;
47const double xmass[5] = {0.000511, 0.105658, 0.139570,0.493677, 0.938272};
48// e mu pi K p
49
50typedef std::vector<int> Vint;
51typedef std::vector<HepLorentzVector> Vp4;
52
53/////////////////////////////////////////////////////////////////////////////
54// //
55// Ks + X //
56// |-->pi+ pi- //
57// //
58// xugf 2009.06 //
59/////////////////////////////////////////////////////////////////////////////
60DECLARE_COMPONENT(inclks)
61inclks::inclks(const std::string& name, ISvcLocator* pSvcLocator) :
62 Algorithm(name, pSvcLocator) {
63
64 m_tuple1 = 0;
65 for(int i = 0; i < 10; i++) m_pass[i] = 0;
66
67//Declare the properties
68 declareProperty("Vr0cut", m_vr0cut=10.0);
69 declareProperty("Vz0cut", m_vz0cut=50.0);
70 declareProperty("CheckDedx", m_checkDedx = 1);
71 declareProperty("CheckTof", m_checkTof = 1);
72
73}
74// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
75StatusCode inclks::initialize(){
76 MsgStream log(msgSvc(), name());
77
78 log << MSG::INFO << "in initialize()" << endmsg;
79
80 StatusCode status;
81
82 if(service("THistSvc", m_thsvc).isFailure()) {
83 log << MSG::ERROR << "Couldn't get THistSvc" << endreq;
84 return StatusCode::FAILURE;
85 }
86 // "DQAHist" is fixed
87 TH1F* inclks_mass = new TH1F("InclKs_mass","INCLUSIVE_Ks_MASS",70,0.46,0.53);
88 inclks_mass->GetXaxis()->SetTitle("M_{#pi#pi} (GeV)");
89 inclks_mass->GetYaxis()->SetTitle("Nentries/0.1MeV/c^{2}");
90 if(m_thsvc->regHist("/DQAHist/InclKs/InclKs_mass", inclks_mass).isFailure()) {
91 log << MSG::ERROR << "Couldn't register inclks_mass" << endreq;
92 }
93
94//*****************************************
95
96 NTuplePtr nt1(ntupleSvc(), "DQAFILE/InclKs");
97 if ( nt1 ) m_tuple1 = nt1;
98 else {
99 m_tuple1 = ntupleSvc()->book ("DQAFILE/InclKs", CLID_ColumnWiseTuple, "ksklx Ntuple");
100 if ( m_tuple1 ) {
101 status = m_tuple1->addItem ("npip", m_npip);
102 status = m_tuple1->addItem ("npim", m_npim);
103 status = m_tuple1->addItem ("ctau", m_ctau);
104 status = m_tuple1->addItem ("lxyz", m_lxyz);
105 status = m_tuple1->addItem ("elxyz", m_elxyz);
106 status = m_tuple1->addItem ("kschi", m_chis);
107 status = m_tuple1->addItem ("mksraw", m_ksRawMass);
108 status = m_tuple1->addItem ("pk0", m_pk0);
109
110 }
111 else {
112 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple1) << endmsg;
113 return StatusCode::FAILURE;
114 }
115 }
116
117 //
118 //--------end of book--------
119 //
120
121 log << MSG::INFO << "successfully return from initialize()" <<endmsg;
122 return StatusCode::SUCCESS;
123
124}
125
126// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
127StatusCode inclks::execute() {
128 StatusCode sc = StatusCode::SUCCESS;
129
130 MsgStream log(msgSvc(), name());
131 log << MSG::INFO << "in execute()" << endreq;
132
133 // DQA
134 // Add the line below at the beginning of execute()
135 //
136 setFilterPassed(false);
137
138 m_pass[0] += 1;
139
140 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
141 SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc(), EventModel::EvtRec::EvtRecEvent);
142 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol(eventSvc(), EventModel::EvtRec::EvtRecTrackCol);
143
144 Vint iks, ipip, ipim, iGood;
145 iGood.clear();
146 iks.clear();
147 ipip.clear();
148 ipim.clear();
149
150 Vp4 ppip, ppim;
151 ppip.clear();
152 ppim.clear();
153
154 int TotCharge = 0;
155 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
156 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
157 if(!(*itTrk)->isMdcTrackValid()) continue;
158 RecMdcTrack *mdcTrk = (*itTrk)->mdcTrack();
159 if(fabs(mdcTrk->z()) >= m_vz0cut) continue;
160 if(mdcTrk->r() >= m_vr0cut) continue;
161 iGood.push_back(i);
162 TotCharge += mdcTrk->charge();
163 }
164 //
165 // Finish Good Charged Track Selection
166 //
167 int nGood = iGood.size();
168
169 //
170 // Charge track number cut
171 //
172
173 if((nGood < 2) || (TotCharge!=0)) return sc;
174
175 m_pass[1] += 1;
176 //
177 // Assign 4-momentum to each charged track
178 //
180 for(int i = 0; i < nGood; i++) {
181 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
182 // if(pid) delete pid;
183 pid->init();
184 pid->setMethod(pid->methodProbability());
185 pid->setChiMinCut(4);
186 pid->setRecTrack(*itTrk);
187 pid->usePidSys(pid->useDedx() | pid->useTof1() | pid->useTof2()); // use PID sub-system
188 pid->identify(pid->onlyPion() | pid->onlyKaon()); // seperater Pion/Kaon
189// pid->identify(pid->onlyPion());
190// pid->identify(pid->onlyKaon());
191 pid->calculate();
192 if(!(pid->IsPidInfoValid())) continue;
193 if(!(*itTrk)->isMdcKalTrackValid()) continue ;
194 RecMdcKalTrack* mdcKalTrk = (*itTrk)->mdcKalTrack();
195// RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
196
197 if(pid->probPion() < 0.001 || (pid->probPion() < pid->probKaon())) continue;
198// if(pid->probPion() < 0.001) continue;
200 HepLorentzVector ptrk;
201 ptrk.setPx(mdcKalTrk->px());
202 ptrk.setPy(mdcKalTrk->py());
203 ptrk.setPz(mdcKalTrk->pz());
204 double p3 = ptrk.mag();
205 ptrk.setE(sqrt(p3*p3+xmass[2]*xmass[2]));
206
207 if(mdcKalTrk->charge() >0) {
208 ipip.push_back(iGood[i]);
209 ppip.push_back(ptrk);
210 } else {
211 ipim.push_back(iGood[i]);
212 ppim.push_back(ptrk);
213 }
214 }
215
216 m_pass[2] += 1;
217 int npip = ipip.size();
218 int npim = ipim.size();
219 m_npip=npip;
220 m_npim=npim;
221
222 if(npip < 1 || npim <1) return sc;
223
224 m_pass[3] += 1;
225
226//
227//****** Second Vertex Check************
228//
229 double chi, delm;
230 double chisq=999.;
231// double delm=100.;
232
233 HepPoint3D vx(0., 0., 0.);
234 HepSymMatrix Evx(3, 0);
235 double bx = 1E+6;
236 double by = 1E+6;
237 double bz = 1E+6;
238 Evx[0][0] = bx*bx;
239 Evx[1][1] = by*by;
240 Evx[2][2] = bz*bz;
241 VertexParameter vxpar;
242 vxpar.setVx(vx);
243 vxpar.setEvx(Evx);
244
245// HepPoint3D bvx(0., 0., 0.);
246// HepSymMatrix bEvx(3, 0);
247// bEvx[0][0] = 0.038*0.038;
248// bEvx[1][1] = 0.00057*0.00057;
249// bEvx[2][2] = 1.5*1.5;
250// VertexParameter bs;
251// bs.setVx(bvx);
252// bs.setEvx(bEvx);
253
254 int ip1=-1, ip2=-1;
255
256 VertexFit *vtxfit0 = VertexFit::instance();
258
259 for(int i1 = 0; i1 < ipip.size(); i1++) {
260 RecMdcKalTrack *pipTrk = (*(evtRecTrkCol->begin()+ipip[i1]))->mdcKalTrack();
262 WTrackParameter wpiptrk(xmass[2], pipTrk->helix(), pipTrk->err());
263
264 for(int i2 = 0; i2 < ipim.size(); i2++) {
265 RecMdcKalTrack *pimTrk = (*(evtRecTrkCol->begin()+ipim[i2]))->mdcKalTrack();
267 WTrackParameter wpimtrk(xmass[2], pimTrk->helix(), pimTrk->err());
268
269 vtxfit0->init();
270 vtxfit0->AddTrack(0, wpiptrk);
271 vtxfit0->AddTrack(1, wpimtrk);
272 vtxfit0->AddVertex(0, vxpar, 0, 1);
273 if(!(vtxfit0->Fit(0))) continue;
274 vtxfit0->BuildVirtualParticle(0);
275
276 vtxfit->init();
277 vtxfit->AddTrack(0, vtxfit0->wVirtualTrack(0));
278 vtxfit->setVpar(vtxfit0->vpar(0));
279 if(!vtxfit->Fit()) continue;
280 chi = vtxfit->chisq();
281
282// if(fabs((vtxfit->p4par().m())-mk0) > delm) continue;
283 if(chi > chisq) continue;
284 delm = fabs((vtxfit->p4par().m())-mk0);
285 chisq = chi;
286 ip1=ipip[i1];
287 ip2=ipim[i2];
288 }
289 }
290//
291// Kshort check
292//
293 if(ip1==-1 || ip2==-1) return sc;
294 m_pass[4] += 1;
295
296 RecMdcKalTrack *pi1Trk = (*(evtRecTrkCol->begin()+ip1))->mdcKalTrack();
298 WTrackParameter wpi1trk(xmass[2], pi1Trk->helix(), pi1Trk->err());
299
300 RecMdcKalTrack *pi2Trk = (*(evtRecTrkCol->begin()+ip2))->mdcKalTrack();
302 WTrackParameter wpi2trk(xmass[2], pi2Trk->helix(), pi2Trk->err());
303 vtxfit0->init();
304 vtxfit0->AddTrack(0, wpi1trk);
305 vtxfit0->AddTrack(1, wpi2trk);
306 vtxfit0->AddVertex(0, vxpar, 0, 1);
307 if(!(vtxfit0->Fit(0))) return sc;
308 vtxfit0->BuildVirtualParticle(0);
309
310 vtxfit->init();
311 vtxfit->AddTrack(0, vtxfit0->wVirtualTrack(0));
312 vtxfit->setVpar(vtxfit0->vpar(0));
313 if(!vtxfit->Fit()) return sc;
314
315 m_ksRawMass = vtxfit->p4par().m();
316 m_ctau = vtxfit->ctau();
317 m_lxyz = vtxfit->decayLength();
318 m_elxyz = vtxfit->decayLengthError();
319 m_chis = vtxfit->chisq();
320 m_pk0 = vtxfit->p4par().rho();
321
322 // DQA
323 if(m_lxyz>0.4 && m_chis<20.){
324
325 TH1 *h(0);
326 if (m_thsvc->getHist("/DQAHist/InclKs/InclKs_mass", h).isSuccess()) {
327 h->Fill(m_ksRawMass);
328 } else {
329 log << MSG::ERROR << "Couldn't retrieve inclks_mass" << endreq;
330 }
331 }
332
333 m_tuple1->write();
334////////////////////////////////////
335//DQA
336// set tag and quality
337 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
338// (*(evtRecTrkCol->begin()+ip1))->setPartId(3);
339// (*(evtRecTrkCol->begin()+ip2))->setPartId(3);
340 (*(evtRecTrkCol->begin()+ip1))->tagPion();
341 (*(evtRecTrkCol->begin()+ip2))->tagPion();
342 // Quality: defined by whether dE/dx or TOF is used to identify particle
343 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
344 // 1 - only dE/dx (can be used for TOF calibration)
345 // 2 - only TOF (can be used for dE/dx calibration)
346 // 3 - Both dE/dx and TOF
347 (*(evtRecTrkCol->begin()+ip1))->setQuality(2);
348 (*(evtRecTrkCol->begin()+ip2))->setQuality(2);
349//--------------------------------------------------
350 // Add the line below at the end of execute(), (before return)
351 //
352 setFilterPassed(true);
353
354 return StatusCode::SUCCESS;
355}
356
357// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
358StatusCode inclks::finalize() {
359
360 MsgStream log(msgSvc(), name());
361 log << MSG::INFO << "in finalize()" << endmsg;
362 log << MSG::INFO << "Total Entries : " << m_pass[0] << endreq;
363 log << MSG::INFO << "Qtot Cut : " << m_pass[1] << endreq;
364 log << MSG::INFO << "PID : " << m_pass[2] << endreq;
365 log << MSG::INFO << "NPI : " << m_pass[3] << endreq;
366 log << MSG::INFO << "Second Vertex Cut : " << m_pass[4] << endreq;
367 return StatusCode::SUCCESS;
368}
369
EvtRecTrackCol::iterator EvtRecTrackIterator
Definition: EvtRecTrack.h:111
std::vector< HepLorentzVector > Vp4
Definition: Gam4pikp.cxx:53
const double xmass[5]
Definition: Gam4pikp.cxx:50
std::vector< int > Vint
Definition: Gam4pikp.cxx:52
NTuple::Tuple * m_tuple1
Definition: MdcHistItem.h:45
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
const HepVector & helix() const
const double px() const
static void setPidType(PidType pidType)
const double pz() const
const double py() const
const HepSymMatrix & err() const
const int charge() const
const double r() const
Definition: DstMdcTrack.h:64
const int charge() const
Definition: DstMdcTrack.h:53
const double z() const
Definition: DstMdcTrack.h:63
int useTof2() const
int methodProbability() const
int useDedx() const
int onlyKaon() const
int onlyPion() const
int useTof1() const
void setChiMinCut(const double chi=4)
void setRecTrack(EvtRecTrack *trk)
double probKaon() const
Definition: ParticleID.h:124
void setMethod(const int method)
Definition: ParticleID.h:94
void identify(const int pidcase)
Definition: ParticleID.h:103
void usePidSys(const int pidsys)
Definition: ParticleID.h:97
static ParticleID * instance()
Definition: ParticleID.cxx:22
bool IsPidInfoValid() const
double probPion() const
Definition: ParticleID.h:123
void calculate()
Definition: ParticleID.cxx:97
void init()
Definition: ParticleID.cxx:27
HepLorentzVector p4par() const
double ctau() const
double decayLength() const
double decayLengthError() const
static SecondVertexFit * instance()
void setVpar(const VertexParameter vpar)
double chisq() const
void AddTrack(const int number, const double mass, const RecMdcTrack *trk)
Definition: TrackPool.cxx:22
WTrackParameter wVirtualTrack(int n) const
Definition: VertexFit.h:92
void init()
Definition: VertexFit.cxx:29
void AddVertex(int number, VertexParameter vpar, std::vector< int > lis)
Definition: VertexFit.cxx:89
static VertexFit * instance()
Definition: VertexFit.cxx:15
VertexParameter vpar(int n) const
Definition: VertexFit.h:89
void BuildVirtualParticle(int number)
Definition: VertexFit.cxx:619
bool Fit()
Definition: VertexFit.cxx:301
void setEvx(const HepSymMatrix &eVx)
void setVx(const HepPoint3D &vx)
Definition: inclks.h:15
StatusCode finalize()
Definition: inclks.cxx:358
StatusCode initialize()
Definition: inclks.cxx:75
StatusCode execute()
Definition: inclks.cxx:127
HepGeom::Point3D< double > HepPoint3D
Definition: inclks.cxx:32
std::vector< HepLorentzVector > Vp4
Definition: inclks.cxx:51
const double xmass[5]
Definition: inclks.cxx:47
std::vector< int > Vint
Definition: inclks.cxx:50
const double mk0
Definition: inclks.cxx:46
_EXTERN_ std::string EvtRecEvent
Definition: EventModel.h:116
_EXTERN_ std::string EvtRecTrackCol
Definition: EventModel.h:117
float ptrk
double double * p3
Definition: qcdloop1.h:76