BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
inclphi Class Reference

#include <inclphi.h>

Inheritance diagram for inclphi:

Public Member Functions

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

Detailed Description

Definition at line 15 of file inclphi.h.

Constructor & Destructor Documentation

◆ inclphi()

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

Definition at line 60 of file inclphi.cxx.

61 : Algorithm( name, pSvcLocator ) {
62
63 m_tuple1 = 0;
64 m_tuple2 = 0;
65 for ( int i = 0; i < 10; i++ ) m_pass[i] = 0;
66
67 // Declare the properties
68 declareProperty( "Vr0cut", m_vr0cut = 5.0 );
69 declareProperty( "Vz0cut", m_vz0cut = 10.0 );
70 declareProperty( "CheckDedx", m_checkDedx = 1 );
71 declareProperty( "CheckTof", m_checkTof = 1 );
72}

Referenced by inclphi().

Member Function Documentation

◆ execute()

StatusCode inclphi::execute ( )

Definition at line 142 of file inclphi.cxx.

142 {
143 StatusCode sc = StatusCode::SUCCESS;
144
145 MsgStream log( msgSvc(), name() );
146 log << MSG::INFO << "in execute()" << endmsg;
147
148 // DQA
149 // Add the line below at the beginning of execute()
150 //
151 setFilterPassed( false );
152
153 m_pass[0] += 1;
154
155 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
156 int run = eventHeader->runNumber();
157 int event = eventHeader->eventNumber();
158
159 SmartDataPtr<EvtRecEvent> evtRecEvent( eventSvc(), EventModel::EvtRec::EvtRecEvent );
160 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol( eventSvc(), EventModel::EvtRec::EvtRecTrackCol );
161
162 m_pass[1] += 1;
163
164 Vint ikp, ikm, iGood;
165 iGood.clear();
166 ikp.clear();
167 ikm.clear();
168
169 Vp4 pkp, pkm;
170 pkp.clear();
171 pkm.clear();
172
173 int TotCharge = 0;
174 for ( int i = 0; i < evtRecEvent->totalCharged(); i++ )
175 {
176 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
177 if ( !( *itTrk )->isMdcTrackValid() ) continue;
178 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
179 if ( fabs( mdcTrk->z() ) >= m_vz0cut ) continue;
180 if ( mdcTrk->r() >= m_vr0cut ) continue;
181 iGood.push_back( i );
182 TotCharge += mdcTrk->charge();
183 }
184 //
185 // Finish Good Charged Track Selection
186 //
187 int nGood = iGood.size();
188
189 //
190 // Charge track number cut
191 //
192
193 if ( ( nGood < 2 ) || ( TotCharge != 0 ) ) return sc;
194
195 m_pass[2] += 1;
196
197 //
198 // Assign 4-momentum to each charged track
199 //
200 ParticleID* pid = ParticleID::instance();
201 for ( int i = 0; i < nGood; i++ )
202 {
203 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
204 // if(pid) delete pid;
205 pid->init();
206 pid->setMethod( pid->methodProbability() );
207 pid->setChiMinCut( 4 );
208 pid->setRecTrack( *itTrk );
209 pid->usePidSys( pid->useDedx() | pid->useTof1() | pid->useTof2() ); // use PID sub-system
210 pid->identify( pid->onlyPion() | pid->onlyKaon() ); // seperater Pion/Kaon
211 pid->calculate();
212 if ( !( pid->IsPidInfoValid() ) ) continue;
213 if ( !( *itTrk )->isMdcKalTrackValid() ) continue;
214 RecMdcKalTrack* mdcKalTrk = ( *itTrk )->mdcKalTrack();
215 // RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
216
217 if ( pid->probKaon() < 0.001 || ( pid->probKaon() < pid->probPion() ) ) continue;
218 // if(pid->probPion() < 0.001) continue;
219 mdcKalTrk->setPidType( RecMdcKalTrack::kaon );
220 HepLorentzVector ptrk;
221 ptrk.setPx( mdcKalTrk->px() );
222 ptrk.setPy( mdcKalTrk->py() );
223 ptrk.setPz( mdcKalTrk->pz() );
224 double p3 = ptrk.mag();
225 ptrk.setE( sqrt( p3 * p3 + mk * mk ) );
226 if ( mdcKalTrk->charge() > 0 )
227 {
228 ikp.push_back( iGood[i] );
229 pkp.push_back( ptrk );
230 }
231 else
232 {
233 ikm.push_back( iGood[i] );
234 pkm.push_back( ptrk );
235 }
236 }
237
238 m_pass[4] += 1;
239 int nkp = ikp.size();
240 int nkm = ikm.size();
241
242 m_nkp = nkp;
243 m_nkm = nkm;
244 m_ncharge = nGood;
245
246 if ( nkp < 1 || nkm < 1 ) return sc;
247
248 m_pass[5] += 1;
249
250 //
251 //**************** Phi Finding ************
252 //
253 //
254
255 HepLorentzVector pphi, pTot;
256 Vint iphi;
257 iphi.clear();
258
259 double difchi0 = 99999.0;
260 int ixk1 = -1;
261 int ixk2 = -1;
262
263 for ( int i = 0; i < nkm; i++ )
264 {
265 for ( int j = 0; j < nkp; j++ )
266 {
267
268 pphi = pkm[i] + pkp[j];
269 m_mphiall = pphi.m();
270 m_pphiall = pphi.rho();
271 m_tuple1->write();
272
273 double difchi = fabs( pphi.m() - mphi );
274 if ( difchi < difchi0 )
275 {
276 difchi0 = difchi;
277 ixk1 = i;
278 ixk2 = j;
279 } // good phi
280 } // K+
281 } // K-
282
283 m_difchi = difchi0;
284
285 if ( ixk1 == -1 ) return sc;
286
287 m_pass[6] += 1;
288
289 pTot = pkm[ixk1] + pkp[ixk2];
290
291 m_pk1 = pkm[ixk1].m();
292 m_mphi = pTot.m();
293 m_pphi = pTot.rho();
294
295 TH1* h( 0 );
296 if ( m_thsvc->getHist( "/DQAHist/InclPhi/InclPhi_mass", h ).isSuccess() )
297 { h->Fill( pTot.m() ); }
298 else { log << MSG::ERROR << "Couldn't retrieve inclphi_mass" << endmsg; }
299 m_tuple2->write();
300 ////////////////////////////////////
301 // DQA
302 // set tag and quality
303 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
304 // (*(evtRecTrkCol->begin()+ikm[ixk1]))->setPartId(4);
305 // (*(evtRecTrkCol->begin()+ikp[ixk2]))->setPartId(4);
306 ( *( evtRecTrkCol->begin() + ikm[ixk1] ) )->tagKaon();
307 ( *( evtRecTrkCol->begin() + ikp[ixk2] ) )->tagKaon();
308 // Quality: defined by whether dE/dx or TOF is used to identify particle
309 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
310 // 1 - only dE/dx (can be used for TOF calibration)
311 // 2 - only TOF (can be used for dE/dx calibration)
312 // 3 - Both dE/dx and TOF
313 ( *( evtRecTrkCol->begin() + ikm[ixk1] ) )->setQuality( 3 );
314 ( *( evtRecTrkCol->begin() + ikp[ixk2] ) )->setQuality( 3 );
315 //--------------------------------------------------
316 // Add the line below at the end of execute(), (before return)
317 //
318 setFilterPassed( true );
319
320 return StatusCode::SUCCESS;
321}
EvtRecTrackCol::iterator EvtRecTrackIterator
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:38
const double mk
Definition Gam4pikp.cxx:33
std::vector< int > Vint
Definition Gam4pikp.cxx:37
IMessageSvc * msgSvc()
static ParticleID * instance()
bool IsPidInfoValid() const
void calculate()
void init()
const double mphi
Definition inclphi.cxx:44

◆ finalize()

StatusCode inclphi::finalize ( )

Definition at line 324 of file inclphi.cxx.

324 {
325
326 MsgStream log( msgSvc(), name() );
327 log << MSG::INFO << "in finalize()" << endmsg;
328 log << MSG::INFO << "Total Entries : " << m_pass[0] << endmsg;
329 log << MSG::INFO << "dstCol, trkCol: " << m_pass[1] << endmsg;
330 log << MSG::INFO << "Ncharge Cut : " << m_pass[2] << endmsg;
331 log << MSG::INFO << "TOF dEdx : " << m_pass[3] << endmsg;
332 log << MSG::INFO << "PID : " << m_pass[4] << endmsg;
333 log << MSG::INFO << "NK Cut : " << m_pass[5] << endmsg;
334 log << MSG::INFO << "Nphi select : " << m_pass[6] << endmsg;
335 return StatusCode::SUCCESS;
336}

◆ initialize()

StatusCode inclphi::initialize ( )

Definition at line 75 of file inclphi.cxx.

75 {
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 {
84 log << MSG::ERROR << "Couldn't get THistSvc" << endmsg;
85 return StatusCode::FAILURE;
86 }
87 // "DQAHist" is fixed
88 TH1F* inclphi_mass = new TH1F( "InclPhi_mass", "INCLUSIVE_PHI_MASS", 80, 1.01, 1.05 );
89 inclphi_mass->GetXaxis()->SetTitle( "M_{KK} (GeV)" );
90 inclphi_mass->GetYaxis()->SetTitle( "Nentries/0.5MeV/c^{2}" );
91 if ( m_thsvc->regHist( "/DQAHist/InclPhi/InclPhi_mass", inclphi_mass ).isFailure() )
92 { log << MSG::ERROR << "Couldn't register inclphi_mass" << endmsg; }
93
94 //*****************************************
95
96 NTuplePtr nt1( ntupleSvc(), "DQAFILE/InclPhi" );
97 if ( nt1 ) m_tuple1 = nt1;
98 else
99 {
100 m_tuple1 = ntupleSvc()->book( "DQAFILE/InclPhi1", CLID_ColumnWiseTuple, "inclphi Ntuple" );
101 if ( m_tuple1 )
102 {
103 status = m_tuple1->addItem( "mphiall", m_mphiall );
104 status = m_tuple1->addItem( "mpphiall", m_pphiall );
105 }
106 else
107 {
108 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple1 ) << endmsg;
109 return StatusCode::FAILURE;
110 }
111 }
112 NTuplePtr nt2( ntupleSvc(), "DQAFILE/InclPhi2" );
113 if ( nt2 ) m_tuple2 = nt2;
114 else
115 {
116 m_tuple2 = ntupleSvc()->book( "DQAFILE/InclPhi2", CLID_ColumnWiseTuple, "inclphi Ntuple" );
117 if ( m_tuple2 )
118 {
119 status = m_tuple2->addItem( "nkp", m_nkp );
120 status = m_tuple2->addItem( "nkm", m_nkm );
121 status = m_tuple2->addItem( "ncharge", m_ncharge );
122 status = m_tuple2->addItem( "difchi0", m_difchi );
123 status = m_tuple2->addItem( "mmphi", m_mphi );
124 status = m_tuple2->addItem( "pphi", m_pphi );
125 status = m_tuple2->addItem( "pk1", m_pk1 );
126 }
127 else
128 {
129 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple2 ) << endmsg;
130 return StatusCode::FAILURE;
131 }
132 }
133 //
134 //--------end of book--------
135 //
136
137 log << MSG::INFO << "successfully return from initialize()" << endmsg;
138 return StatusCode::SUCCESS;
139}
INTupleSvc * ntupleSvc()

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