StFemtoDstInclusiveSelector.cxx 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593
  1. #include "StFemtoDstInclusiveSelector.h"
  2. #include "StBTofHeader.h"
  3. #include "TBranch.h"
  4. #include "StPhysicalHelixD.hh"
  5. ClassImp(StFemtoDstInclusiveSelector)
  6. // Set maximum file size to 1.9 GB (Root has a 2GB limit)
  7. #define MAXFILESIZE 1900000000
  8. //_________________
  9. StFemtoDstInclusiveSelector::StFemtoDstInclusiveSelector(StMuDstMaker *muMaker,
  10. const Char_t *oFileName) {
  11. std::cout << "StFemtoDstInclusiveSelector::Constructor - Creating an instance...";
  12. mOutFileName = oFileName;
  13. mMuDstMaker = muMaker;
  14. mMuDst = mMuDstMaker->muDst();
  15. mRefMultCorrUtil = NULL;
  16. mCompression = 9;
  17. mCollisionTypeAuAu = false; // pp: false; AuAu: true
  18. mLoopAllVert = false; // Not first primary vertex does not contain tracks with fast detectors (TOF)
  19. mEventIsGood = false;
  20. mNEventsIn = 0;
  21. mNEventsPassed = 0;
  22. mAuAuCorrectZdc = true;
  23. mIsGoodTrack = false;
  24. mCurrentTrigger = 0;
  25. //Initialize event cut variables
  26. mPrimVtxZ[0] = -70.;
  27. mPrimVtxZ[1] = 70.;
  28. mPrimVtxR[0] = -1.;
  29. mPrimVtxR[1] = 10.;
  30. mPrimVtxVpdVzDiff[0] = -5.;
  31. mPrimVtxVpdVzDiff[1]= 5.;
  32. mPrimVtxXShift = 0.;
  33. mPrimVtxYShift = 0.;
  34. //Initialize single-particle cut variables
  35. mTrackMomentum[0] = 0.1;
  36. mTrackMomentum[1] = 1.7;
  37. mTrackDca[0] = 0.;
  38. mTrackDca[1] = 5.;
  39. mTrackDcaGlobal[0] = 0.;
  40. mTrackDcaGlobal[1] = 5.;
  41. mTrackNHits[0] = 14;
  42. mTrackNHits[1] = 50;
  43. mTrackEta[0] = -1.;
  44. mTrackEta[1] = 1.;
  45. mTrackFlag[0] = 0;
  46. mTrackFlag[1] = 1000;
  47. mSelectPions = true;
  48. mSelectTpcPions = true;
  49. mPionTpcMomentum[0] = 0.15; mPionTpcMomentum[1] = 0.65;
  50. mPionTpcNSigmaElectron[0]= -2.; mPionTpcNSigmaElectron[1]= 2.;
  51. mPionTpcNSigmaPion[0]= -3.; mPionTpcNSigmaPion[1]= 3.;
  52. mPionTpcNSigmaKaon[0]= -1.5; mPionTpcNSigmaKaon[1]= 1.5;
  53. mPionTpcNSigmaProton[0]= -3.; mPionTpcNSigmaProton[1]= 3.;
  54. mSelectTofPions = true;
  55. mPionTofMomentum[0] = 0.15; mPionTofMomentum[1] = 1.75;
  56. mPionMassSqr[0] = -0.01; mPionMassSqr[1] = 0.09;
  57. mSelectKaons = true;
  58. mSelectTpcKaons = true;
  59. mKaonTpcMomentum[0] = 0.15; mKaonTpcMomentum[1] = 0.65;
  60. mKaonTpcNSigmaElectron[0]= -2.; mKaonTpcNSigmaElectron[1]= 2.;
  61. mKaonTpcNSigmaPion[0]= -1.5; mKaonTpcNSigmaPion[1]= 1.5;
  62. mKaonTpcNSigmaKaon[0]= -3.; mKaonTpcNSigmaKaon[1]= 3.;
  63. mKaonTpcNSigmaProton[0]= -3.; mKaonTpcNSigmaProton[1]= 3.;
  64. mSelectTofKaons = true;
  65. mKaonTofMomentum[0] = 0.15; mKaonTofMomentum[1] = 1.75;
  66. mKaonMassSqr[0] = 0.16; mKaonMassSqr[1] = 0.35;
  67. mSelectProtons = false;
  68. mSelectTpcProtons = true;
  69. mProtonTpcMomentum[0] = 0.15; mProtonTpcMomentum[1] = 1.2;
  70. mProtonTpcNSigmaElectron[0]= -3.; mProtonTpcNSigmaElectron[1]= 3.;
  71. mProtonTpcNSigmaPion[0]= -3.; mProtonTpcNSigmaPion[1]= 3.;
  72. mProtonTpcNSigmaKaon[0]= -3.; mProtonTpcNSigmaKaon[1]= 3.;
  73. mProtonTpcNSigmaProton[0]= -3.; mProtonTpcNSigmaProton[1]= 3.;
  74. mSelectTofProtons = true;
  75. mProtonTofMomentum[0] = 0.15; mProtonTofMomentum[1] = 1.75;
  76. mProtonMassSqr[0] = 0.65; mProtonMassSqr[1] = 1.05;
  77. std::cout << "\t[DONE]" << std::endl;
  78. }
  79. //_________________
  80. StFemtoDstInclusiveSelector::~StFemtoDstInclusiveSelector() {
  81. /* do nothing */
  82. }
  83. //_________________
  84. Int_t StFemtoDstInclusiveSelector::Init() {
  85. std::cout << "StFemtoDstInclusiveSelector::Init - Initializing the maker"
  86. << std::endl
  87. << "Creating output file: " << mOutFileName;
  88. mFemtoEvent = new StFemtoEvent();
  89. mOutFile = new TFile(mOutFileName, "RECREATE");
  90. mOutFile->SetCompressionLevel(mCompression);
  91. std::cout << "\t[DONE]" << std::endl;
  92. mTree = new TTree("StFemtoDst","StFemtoDst");
  93. mTree->SetMaxTreeSize(MAXFILESIZE);
  94. //mTree->SetAutoSave(1000000);
  95. mTree->Branch("StFemtoEvent","StFemtoEvent", &mFemtoEvent);
  96. if(mCollisionTypeAuAu) { //For AuAu only
  97. if(!mRefMultCorrUtil) {
  98. std::cout << "StFemtoDstInclusiveSelector::Init - Initializing StRefMultCorr...";
  99. mRefMultCorrUtil = new StRefMultCorr("refmult");
  100. std::cout << "\t[DONE]" << endl;
  101. }
  102. }
  103. mNBytes = 0;
  104. std::cout << "StFemtoDstInclusiveSelector::Init - Initialization has been finished"
  105. << std::endl;
  106. return StMaker::Init();
  107. }
  108. //_________________
  109. void StFemtoDstInclusiveSelector::Clear(Option_t *option) {
  110. StMaker::Clear(option);
  111. }
  112. //_________________
  113. Int_t StFemtoDstInclusiveSelector::Make() {
  114. // Get MuEvent from the MuDst
  115. mNEventsIn++;
  116. StMuEvent *mMuEvent = mMuDst->event();
  117. if (!AcceptTrigger(mMuEvent)) { // If trigger the trigger was not found
  118. return kStOk;
  119. }
  120. if (mMuEvent->refMult() < 0) { // Multiplicity can't be negative
  121. return kStOk;
  122. }
  123. UShort_t mNVertices = mLoopAllVert ? mMuDst->numberOfPrimaryVertices() : 1;
  124. Short_t mNPrimTracks = mMuDst->numberOfPrimaryTracks();
  125. Short_t mNGlobTracks = mMuDst->numberOfGlobalTracks();
  126. //Reasonable amount of tracks
  127. if (mNPrimTracks < 0 || mNPrimTracks > 10000) {
  128. return kStOk;
  129. }
  130. Bool_t mEventIsGood = false;
  131. Float_t mVpdVz = 0.;
  132. StBTofHeader *mTofHeader = mMuDst->btofHeader();
  133. // Centrality variables
  134. Int_t mCent16;
  135. Double_t mCentWeight, mRefMultCorrVal;
  136. // If pp collisions
  137. if (!mCollisionTypeAuAu) {
  138. mCent16 = -1;
  139. mCentWeight = -999.;
  140. }
  141. Float_t mRanking = 0;
  142. StMuPrimaryVertex *mPrimVertex = NULL;
  143. StMuTrack *mPrimTrack = NULL;
  144. StMuTrack *mGlobTrack = NULL;
  145. Float_t mMassSqr = -999.;
  146. Float_t mBeta = -999.;
  147. Bool_t mIsTofTrack = false;
  148. Bool_t mIsTpcPion = false;
  149. Bool_t mIsTofPion = false;
  150. Bool_t mIsTpcKaon = false;
  151. Bool_t mIsTofKaon = false;
  152. Bool_t mIsTpcProton = false;
  153. Bool_t mIsTofProton = false;
  154. //Calculate the amount of the pile-up vertices
  155. UShort_t mNVtxInBunchCross = 0;
  156. for(UShort_t iVert=0; iVert<mMuDst->numberOfPrimaryVertices(); iVert++) {
  157. mPrimVertex = mMuDst->primaryVertex(iVert);
  158. if(!mPrimVertex || !mTofHeader) continue;
  159. if(mPrimVertex->ranking() <= 0.) continue;
  160. if(TMath::Abs(mPrimVertex->position().x()) < 1e-5 &&
  161. TMath::Abs(mPrimVertex->position().y()) < 1e-5 &&
  162. TMath::Abs(mPrimVertex->position().z()) < 1e-5) {
  163. continue;
  164. }
  165. if(TMath::Abs( mPrimVertex->position().z() - mTofHeader->vpdVz()) < 5. ) {
  166. mNVtxInBunchCross++;
  167. }
  168. } //for(unsigned short iVert=0; iVert<mNVertices; iVert++)
  169. //There should be only one primary vertex that caused the trigger
  170. if(mNVtxInBunchCross > 1) {
  171. return kStOk;
  172. }
  173. // Vertex loop
  174. for (UShort_t iVert = 0; iVert < mNVertices; iVert++) {
  175. mEventIsGood = false;
  176. mPrimVertex = mMuDst->primaryVertex(iVert);
  177. if (!mPrimVertex) {
  178. //std::cout << "No primary vertices in the event. Skip..." << std::endl;
  179. continue;
  180. }
  181. mRanking = mPrimVertex->ranking();
  182. //Positive ranking only
  183. if (mPrimVertex->ranking() <= 0) continue;
  184. //Not (0,0,0) position of the primary vertex
  185. if (TMath::Abs(mPrimVertex->position().x()) < 1e-5 &&
  186. TMath::Abs(mPrimVertex->position().y()) < 1e-5 &&
  187. TMath::Abs(mPrimVertex->position().z()) < 1e-5 ) continue;
  188. StThreeVectorF mVertPosition = mPrimVertex->position();
  189. if(mTofHeader) { //In case of the absence of the vpd information
  190. mVpdVz = mTofHeader->vpdVz();
  191. }
  192. else {
  193. mVpdVz = mVertPosition.z();
  194. }
  195. // Event cut
  196. if (!AcceptPrimVtx(mVertPosition, mVpdVz)) continue;
  197. mEventIsGood = true;
  198. // Set StRefMultCorr
  199. if (mCollisionTypeAuAu) { // Only for AuAu collisions
  200. mRefMultCorrUtil->init(mMuEvent->runId());
  201. //Check for the bad run
  202. if (mRefMultCorrUtil->isBadRun(mMuEvent->runId())) {
  203. std::cout << "StRefMultCorrUtil::IsBadRun" << std::endl;
  204. continue;
  205. }
  206. if (mAuAuCorrectZdc == true) {
  207. mRefMultCorrUtil->initEvent(mMuEvent->refMult(),
  208. mVertPosition.z(),
  209. mMuEvent->runInfo().zdcCoincidenceRate());
  210. }
  211. else {
  212. mRefMultCorrUtil->initEvent(mMuEvent->refMult(),
  213. mVertPosition.z());
  214. }
  215. mCent16 = mRefMultCorrUtil->getCentralityBin16();
  216. if (mCent16 < 0) continue; //Some AuAu collisions have -1 value
  217. mCentWeight = mRefMultCorrUtil->getWeight();
  218. mRefMultCorrVal = mRefMultCorrUtil->getRefMultCorr();
  219. }
  220. // Write event information to the FemtoDst
  221. mNEventsPassed++;
  222. mFemtoEvent->SetEventId( mMuEvent->eventId() );
  223. mFemtoEvent->SetRunId( mMuEvent->runId() );
  224. mFemtoEvent->SetCollisionType( mCollisionTypeAuAu );
  225. mFemtoEvent->SetRefMult( mMuEvent->refMult() );
  226. mFemtoEvent->SetRefMultCorr( mRefMultCorrVal );
  227. mFemtoEvent->SetRefMultCorrWeight( mCentWeight );
  228. mFemtoEvent->SetCent16( mCent16 );
  229. mFemtoEvent->SetRefMultPos( mMuEvent->refMultPos() );
  230. mFemtoEvent->SetNumberOfBTofHit( mMuDst->numberOfBTofHit() ); //mMuEvent->btofTrayMultiplicity() );
  231. //mFemtoEvent->SetZDCe( mMuEvent->zdcTriggerDetector().adc(4) );
  232. //mFemtoEvent->SetZDCw( mMuEvent->zdcTriggerDetector().adc(0) );
  233. mFemtoEvent->SetNumberOfPrimaryTracks( mNPrimTracks );
  234. mFemtoEvent->SetNumberOfGlobalTracks( mNGlobTracks );
  235. mFemtoEvent->SetMagneticField( mMuEvent->magneticField() );
  236. mFemtoEvent->SetPrimaryVertexPosition( mVertPosition.x(),
  237. mVertPosition.y(),
  238. mVertPosition.z() );
  239. mFemtoEvent->SetVpdVz(mVpdVz);
  240. mFemtoEvent->SetTriggerId(mCurrentTrigger);
  241. mFemtoEvent->SetPrimaryVertexRanking(mRanking);
  242. UInt_t refMult2 = 0;
  243. UInt_t refMult2Pos = 0;
  244. //Primary track loop
  245. for(UShort_t iTrk = 0; iTrk < mNPrimTracks; iTrk++) {
  246. mPrimTrack = mMuDst->primaryTracks(iTrk);
  247. mGlobTrack = mMuDst->globalTracks(mPrimTrack->index2Global());
  248. //RefMult2 calculation
  249. if ( AcceptRefMult2(mPrimTrack) ) {
  250. refMult2++;
  251. if(mPrimTrack->charge() > 0) refMult2Pos++;
  252. }
  253. //Track cut
  254. if ( !AcceptTrack(mPrimTrack, iVert) ) continue;
  255. mIsTofTrack = false;
  256. mIsTofTrack = IsTofTrack(mGlobTrack);
  257. mIsTpcPion = false;
  258. mIsTofPion = false;
  259. mIsTpcKaon = false;
  260. mIsTofKaon = false;
  261. mIsTpcProton = false;
  262. mIsTofProton = false;
  263. mBeta = -999.;
  264. mMassSqr = -999.;
  265. // Check if track has ToF signal
  266. if(mIsTofTrack) {
  267. mBeta = mGlobTrack->btofPidTraits().beta();
  268. mMassSqr = mPrimTrack->momentum().mag2() * ( 1./(mBeta * mBeta) - 1.);
  269. }
  270. // Pion selection
  271. if (mSelectPions) {
  272. if (mSelectTpcPions)
  273. mIsTpcPion = IsTpcPion(mPrimTrack);
  274. if (mSelectTofPions && mIsTofTrack)
  275. mIsTofPion = IsTofPion(mMassSqr, mGlobTrack->p().mag());
  276. }
  277. // Kaon selection
  278. if (mSelectKaons) {
  279. if (mSelectTpcKaons)
  280. mIsTpcKaon = IsTpcKaon(mPrimTrack);
  281. if (mSelectTofKaons && mIsTofTrack)
  282. mIsTofKaon = IsTofKaon(mMassSqr, mGlobTrack->p().mag());
  283. }
  284. // Proton selection
  285. if (mSelectProtons) {
  286. if (mSelectTpcProtons)
  287. mIsTpcProton = IsTpcProton(mPrimTrack);
  288. if (mSelectTofProtons && mIsTofTrack)
  289. mIsTofProton = IsTofProton(mMassSqr, mGlobTrack->p().mag());
  290. }
  291. // Write StFemtoTrack
  292. if (mIsTpcPion || mIsTofPion || mIsTpcKaon || mIsTofKaon ||
  293. mIsTpcProton || mIsTofProton) {
  294. mFTrack = mFemtoEvent->AddFemtoTrack();
  295. mFTrack->SetId( mPrimTrack->id() );
  296. mFTrack->SetNHits( (Char_t)(mPrimTrack->charge()*mPrimTrack->nHits()) );
  297. mFTrack->SetNHitsPoss( mPrimTrack->nHitsPoss(kTpcId) );
  298. //mFTrack->SetNHitsDedx( mPrimTrack->nHitsDedx() );
  299. mFTrack->SetNSigmaElectron( mPrimTrack->nSigmaElectron() );
  300. mFTrack->SetNSigmaPion( mPrimTrack->nSigmaPion() );
  301. mFTrack->SetNSigmaKaon( mPrimTrack->nSigmaKaon() );
  302. mFTrack->SetNSigmaProton( mPrimTrack->nSigmaProton() );
  303. mFTrack->SetDedx( mPrimTrack->dEdx() );
  304. //mFTrack->SetChi2( mPrimTrack->chi2() );
  305. mFTrack->SetMap0( mGlobTrack->topologyMap().data(0) );
  306. mFTrack->SetMap1( mGlobTrack->topologyMap().data(1) );
  307. mFTrack->SetP( mPrimTrack->momentum().x(),
  308. mPrimTrack->momentum().y(),
  309. mPrimTrack->momentum().z() );
  310. mFTrack->SetDCAxGlobal( mPrimTrack->dcaGlobal().x() );
  311. mFTrack->SetDCAyGlobal( mPrimTrack->dcaGlobal().y() );
  312. mFTrack->SetDCAzGlobal( mPrimTrack->dcaGlobal().z() );
  313. mFTrack->SetGlobalP( mGlobTrack->p().x(),
  314. mGlobTrack->p().y(),
  315. mGlobTrack->p().z() );
  316. if(mIsTofTrack) {
  317. mFTrack->SetBeta( mBeta );
  318. }
  319. else {
  320. mFTrack->SetBeta( -999. );
  321. }
  322. }
  323. } //for(Int_t iTrk = 0; iTrk < mNPrimTracks; iTrk++)
  324. //Store RefMult2
  325. mFemtoEvent->SetRefMult2(refMult2);
  326. mFemtoEvent->SetRefMult2Pos(refMult2Pos);
  327. //
  328. // If event is good then write it
  329. //
  330. if(mEventIsGood) {
  331. mNBytes += mTree->Fill();
  332. mFemtoEvent->Clear();
  333. }
  334. } //(Int_t iVert = 0; iVert < mNVertices; iVert++)
  335. return kStOk;
  336. }
  337. //_________________
  338. Bool_t StFemtoDstInclusiveSelector::AcceptRefMult2(StMuTrack *trk) {
  339. Bool_t mIsGoodTrack = false;
  340. if ( trk->flag()> 0 &&
  341. trk->charge() != 0 &&
  342. trk->nHitsFit() >= 10 &&
  343. trk->dca().mag() < 3.0 &&
  344. trk->momentum().mag() >= 1e-10 &&
  345. TMath::Abs(trk->eta()) <= 1.0) {
  346. mIsGoodTrack = true;
  347. }
  348. return mIsGoodTrack;
  349. }
  350. //_________________
  351. Bool_t StFemtoDstInclusiveSelector::AcceptTrigger(StMuEvent *muEvent) {
  352. //We will store the bit mask for the list of triggers.
  353. //The list of triggers should be the same and in the same
  354. //order, in order to read it later
  355. Bool_t mIsGoodTrigger = false;
  356. mCurrentTrigger = 0;
  357. if(mTriggerIdCollection.empty()) {
  358. mIsGoodTrigger = true;
  359. }
  360. else {
  361. for (UInt_t iTrg = 0; iTrg < mTriggerIdCollection.size(); iTrg++) {
  362. if(muEvent->triggerIdCollection().nominal().isTrigger(mTriggerIdCollection.at(iTrg))) {
  363. mIsGoodTrigger = true;
  364. mCurrentTrigger |= (1<<iTrg);
  365. break;
  366. }
  367. } //for(UInt_t iTrg=0; iTrg<mTriggerIdCollection.size(); iTrg++)
  368. }
  369. return mIsGoodTrigger;
  370. }
  371. //_________________
  372. Bool_t StFemtoDstInclusiveSelector::AcceptPrimVtx(StThreeVectorF vtxPos,
  373. Float_t vpdVz) {
  374. Bool_t mIsGoodVtx = false;
  375. Float_t mVtxX = vtxPos.x() - mPrimVtxXShift;
  376. Float_t mVtxY = vtxPos.y() - mPrimVtxYShift;
  377. Float_t vtxR = TMath::Sqrt(mVtxX*mVtxX + mVtxY*mVtxY);
  378. Float_t vpdDiff = vtxPos.z() - vpdVz;
  379. if( vtxPos.z() >= mPrimVtxZ[0] &&
  380. vtxPos.z() <= mPrimVtxZ[1] &&
  381. vtxR >= mPrimVtxR[0] &&
  382. vtxR <= mPrimVtxR[1] &&
  383. vpdDiff >= mPrimVtxVpdVzDiff[0] &&
  384. vpdDiff <= mPrimVtxVpdVzDiff[1]) {
  385. mIsGoodVtx = true;
  386. }
  387. return mIsGoodVtx;
  388. }
  389. //_________________
  390. Bool_t StFemtoDstInclusiveSelector::AcceptTrack(StMuTrack *trk, UShort_t vtxInd) {
  391. Bool_t mIsGoodTrack = false;
  392. if( trk->vertexIndex() == vtxInd &&
  393. trk->momentum().mag() >= mTrackMomentum[0] && trk->momentum().mag() <= mTrackMomentum[1] &&
  394. trk->type() == 1 &&
  395. trk->eta() >= mTrackEta[0] && trk->eta() <= mTrackEta[1] &&
  396. trk->nHits() >= mTrackNHits[0] && trk->nHits() <= mTrackNHits[1] &&
  397. trk->flag() >= mTrackFlag[0] && trk->flag() <= mTrackFlag[1] &&
  398. trk->dca(vtxInd).perp() >= mTrackDca[0] && trk->dca(vtxInd).perp() <= mTrackDca[1] &&
  399. trk->dcaGlobal(vtxInd).perp() >= mTrackDcaGlobal[0] &&
  400. trk->dcaGlobal(vtxInd).perp() <= mTrackDcaGlobal[1]) {
  401. mIsGoodTrack = true;
  402. }
  403. return mIsGoodTrack;
  404. }
  405. //_________________
  406. Bool_t StFemtoDstInclusiveSelector::IsTofTrack(StMuTrack *trk) { //Only for globals
  407. Bool_t mIsTofHit = false;
  408. if(trk->btofPidTraits().beta() > 0 &&
  409. trk->btofPidTraits().timeOfFlight() > 0) {
  410. mIsTofHit = true;
  411. }
  412. return mIsTofHit;
  413. }
  414. //_________________
  415. Int_t StFemtoDstInclusiveSelector::Finish() {
  416. mOutFile->cd();
  417. mOutFile->Write();
  418. //mTree->Write();
  419. mOutFile->Close();
  420. delete mTree;
  421. delete mOutFile;
  422. std::cout << "*********************************************" << std::endl
  423. << "StFemtoDstInclusiveSelector has been finished" << std::endl
  424. << "\t nEventsPassed : " << mNEventsPassed << std::endl
  425. << "\t nEventsProcessed: " << mNEventsIn << std::endl
  426. << "**********************************************" << std::endl;
  427. return kStOk;
  428. }
  429. //_________________
  430. Bool_t StFemtoDstInclusiveSelector::IsTpcPion(StMuTrack *pTrk) {
  431. Float_t p = pTrk->p().mag();
  432. Float_t nSElectron = pTrk->nSigmaElectron();
  433. Float_t nSPion = pTrk->nSigmaPion();
  434. Float_t nSKaon = pTrk->nSigmaKaon();
  435. Float_t nSProton = pTrk->nSigmaProton();
  436. return ( (p >= mPionTpcMomentum[0]) && (p <= mPionTpcMomentum[1]) &&
  437. ( (nSPion >= mPionTpcNSigmaPion[0]) &&
  438. (nSPion <= mPionTpcNSigmaPion[1]) ) &&
  439. ( (nSElectron <= mPionTpcNSigmaElectron[0]) ||
  440. (nSElectron >= mPionTpcNSigmaElectron[1])) &&
  441. ( (nSKaon <= mPionTpcNSigmaKaon[0]) ||
  442. (nSKaon >= mPionTpcNSigmaKaon[1]) ) &&
  443. ( (nSProton <= mPionTpcNSigmaProton[0]) ||
  444. (nSProton >= mPionTpcNSigmaProton[1]) ) );
  445. }
  446. //_________________
  447. Bool_t StFemtoDstInclusiveSelector::IsTpcKaon(StMuTrack *pTrk) {
  448. Float_t p = pTrk->p().mag();
  449. Float_t nSElectron = pTrk->nSigmaElectron();
  450. Float_t nSPion = pTrk->nSigmaPion();
  451. Float_t nSKaon = pTrk->nSigmaKaon();
  452. Float_t nSProton = pTrk->nSigmaProton();
  453. return ( (p >= mKaonTpcMomentum[0]) && (p <= mKaonTpcMomentum[1]) &&
  454. ( (nSKaon >= mKaonTpcNSigmaKaon[0]) &&
  455. (nSKaon <= mKaonTpcNSigmaKaon[1]) ) &&
  456. ( (nSElectron <= mKaonTpcNSigmaElectron[0]) ||
  457. (nSElectron >= mKaonTpcNSigmaElectron[1]) ) &&
  458. ( (nSPion <= mKaonTpcNSigmaPion[0]) ||
  459. (nSPion >= mKaonTpcNSigmaPion[1]) ) &&
  460. ( (nSProton <= mKaonTpcNSigmaProton[0]) ||
  461. (nSProton >= mKaonTpcNSigmaProton[1]) ) );
  462. }
  463. //_________________
  464. Bool_t StFemtoDstInclusiveSelector::IsTpcProton(StMuTrack *pTrk) {
  465. Float_t p = pTrk->p().mag();
  466. Float_t nSElectron = pTrk->nSigmaElectron();
  467. Float_t nSPion = pTrk->nSigmaPion();
  468. Float_t nSKaon = pTrk->nSigmaKaon();
  469. Float_t nSProton = pTrk->nSigmaProton();
  470. return ( (p >= mProtonTpcMomentum[0]) && (p <= mProtonTpcMomentum[1]) &&
  471. ( (nSProton >= mProtonTpcNSigmaProton[0]) &&
  472. (nSProton <= mProtonTpcNSigmaProton[1]) ) &&
  473. ( (nSElectron <= mProtonTpcNSigmaElectron[0]) ||
  474. (nSElectron >= mProtonTpcNSigmaElectron[1]) ) &&
  475. ( (nSPion <= mProtonTpcNSigmaPion[0]) ||
  476. (nSPion >= mProtonTpcNSigmaPion[1]) ) &&
  477. ( (nSKaon <= mProtonTpcNSigmaKaon[0]) ||
  478. (nSKaon >= mProtonTpcNSigmaKaon[1]) ) );
  479. }
  480. //_________________
  481. Bool_t StFemtoDstInclusiveSelector::IsTofPion(Float_t mSqr, Float_t p) {
  482. return ( (p >= mPionTofMomentum[0]) && (p <= mPionTofMomentum[1]) &&
  483. (mSqr >= mPionMassSqr[0]) && (mSqr <= mPionMassSqr[1]) );
  484. }
  485. //_________________
  486. Bool_t StFemtoDstInclusiveSelector::IsTofKaon(Float_t mSqr, Float_t p) {
  487. return ( (p >= mKaonTofMomentum[0]) && (p <= mKaonTofMomentum[1]) &&
  488. (mSqr >= mKaonMassSqr[0]) && (mSqr <= mKaonMassSqr[1]) );
  489. }
  490. //_________________
  491. Bool_t StFemtoDstInclusiveSelector::IsTofProton(Float_t mSqr, Float_t p) {
  492. return ( (p >= mProtonTofMomentum[0]) && (p <= mProtonTofMomentum[1]) &&
  493. (mSqr >= mProtonMassSqr[0]) && (mSqr <= mProtonMassSqr[1]) );
  494. }
  495. //_________________
  496. void StFemtoDstInclusiveSelector::AddTriggerId(const UInt_t& id) {
  497. mTriggerIdCollection.push_back(id);
  498. }