StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFwdAnalysisMaker.cxx
1 #include "TVector3.h"
2 #include "TH1F.h"
3 #include "TH2F.h"
4 
5 #include "St_base/StMessMgr.h"
6 
7 #include "StEvent/StEvent.h"
9 #include "StEvent/StFcsCluster.h"
10 #include "StEvent/StFttCollection.h"
11 #include "StEvent/StFcsCollection.h"
12 #include "StEvent/StFwdTrack.h"
13 #include "StEvent/StFwdTrackCollection.h"
14 
15 #include "StMuDSTMaker/COMMON/StMuDstMaker.h"
16 #include "StMuDSTMaker/COMMON/StMuDst.h"
17 #include "StMuDSTMaker/COMMON/StMuFwdTrack.h"
18 #include "StMuDSTMaker/COMMON/StMuFwdTrackCollection.h"
19 
20 #include "StFcsDbMaker/StFcsDb.h"
21 #include "StFwdTrackMaker/Common.h"
22 #include "StFwdUtils/StFwdAnalysisMaker.h"
23 
24 //________________________________________________________________________
25 StFwdAnalysisMaker::StFwdAnalysisMaker() : StMaker("fwdAna"){
26  setLocalOutputFile( "" ); // default off
27 };
29 
30  if ( mLocalOutputFile != "" ){
31  auto prevDir = gDirectory;
32 
33  // output file name
34  TFile *fOutput = new TFile(mLocalOutputFile, "RECREATE");
35  fOutput->cd();
36  for (auto nh : mHists) {
37  nh.second->SetDirectory(gDirectory);
38  nh.second->Write();
39  }
40 
41  // restore previous directory
42  gDirectory = prevDir;
43 
44  LOG_INFO << "Done writing StFwdAnalysisMaker output to local file : " << mLocalOutputFile << endm;
45  }
46 
47  return kStOk;
48 }
49 
50 //________________________________________________________________________
51 int StFwdAnalysisMaker::Init() {
52  LOG_DEBUG << "StFwdAnalysisMaker::Init" << endm;
53 
54  AddHist( mHists["fwdMultFailed"] = new TH1F("fwdMultFailed", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
55  AddHist( mHists["fwdMultAll"] = new TH1F("fwdMultAll", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
56  AddHist( mHists["fwdMultGood"] = new TH1F("fwdMultGood", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
57  AddHist( mHists["fwdMultFST"] = new TH1F("fwdMultFST", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
58  AddHist( mHists["nHitsFit"] = new TH1F("nHitsFit", ";nHitsFit; counts", 10, 0, 10) );
59  AddHist( mHists["fwdMultEcalMatch"] = new TH1F("fwdMultEcalMatch", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
60  AddHist( mHists["fwdMultHcalMatch"] = new TH1F("fwdMultHcalMatch", ";N_{ch}^{FWD}; counts", 100, 0, 100) );
61  AddHist( mHists["fwdMultEcalClusters"] = new TH1F("fwdMultEcalClusters", ";N_{Clu}^{ECAL}; counts", 100, 0, 100) );
62  AddHist( mHists["fwdMultHcalClusters"] = new TH1F("fwdMultHcalClusters", ";N_{Clu}^{HCAL}; counts", 100, 0, 100) );
63  AddHist( mHists["eta"] = new TH1F("eta", ";#eta; counts", 100, 0, 5) );
64  AddHist( mHists["phi"] = new TH1F("phi", ";#phi; counts", 100, -3.1415926, 3.1415926) );
65  AddHist( mHists["pt"] = new TH1F("pt", "; pT; counts", 500, 0, 10) );
66  AddHist( mHists["charge"] = new TH1F("charge", "; charge; counts", 4, -2, 2) );
67  AddHist( mHists["ecalMatchPerTrack"] = new TH1F("ecalMatchPerTrack", ";N_{match} / track; counts", 5, 0, 5) );
68  AddHist( mHists["hcalMatchPerTrack"] = new TH1F("hcalMatchPerTrack", ";N_{match} / track; counts", 5, 0, 5) );
69  AddHist( mHists["matchedEcalEnergy"] = new TH1F("matchedEcalEnergy", ";Energy; counts", 100, 0, 15) );
70  AddHist( mHists["matchedHcalEnergy"] = new TH1F("matchedHcalEnergy", ";Energy; counts", 100, 0, 15) );
71  AddHist( mHists["ecalEnergy"] = new TH1F("ecalEnergy", ";Energy; counts", 100, 0, 15) );
72  AddHist( mHists["hcalEnergy"] = new TH1F("hcalEnergy", ";Energy; counts", 100, 0, 15) );
73  AddHist( mHists["ecalXY"] = new TH2F( "ecalXY", ";ecalX;ecalY", 200, -200, 200, 200, -200, 200 ) );
74  AddHist( mHists["hcalXY"] = new TH2F( "hcalXY", ";hcalX;hcalY", 200, 0, 50, 200, 0, 50 ) );
75  AddHist( mHists["ecaldX"] = new TH1F( "ecaldX", ";dx (trk - ecal); counts", 400, -200, 200 ) );
76  AddHist( mHists["matchedEcaldX"] = new TH1F( "matchedEcaldX", ";dx (trk - ecal); counts", 400, -200, 200 ) );
77  AddHist( mHists["ecaldY"] = new TH1F( "ecaldY", ";dy (trk - ecal); counts", 400, -200, 200 ) );
78  AddHist( mHists["matchedEcaldY"] = new TH1F( "matchedEcaldY", ";dy (trk - ecal); counts", 400, -200, 200 ) );
79  AddHist( mHists["ecaldR"] = new TH1F( "ecaldR", ";dr (trk - ecal); counts", 400, 0, 400 ) );
80  AddHist( mHists["ecalMindR"] = new TH1F( "ecalMindR", ";dr (trk - ecal); counts", 400, 0, 400 ) );
81  AddHist( mHists["matchedEcaldR"] = new TH1F( "matchedEcaldR", ";dr (trk - ecal); counts", 400, 0, 400 ) );
82  AddHist( mHists["hcaldX"] = new TH1F( "hcaldX", ";dx (trk - hcal); counts", 400, -200, 200 ) );
83  AddHist( mHists["hcaldXdNFit"] = new TH2F( "hcaldXdNFit", ";dx (trk - hcal); nFit", 400, -200, 200, 10, 0, 10 ) );
84  AddHist( mHists["matchedHcaldX"] = new TH1F( "matchedHcaldX", ";dx (trk - hcal); counts", 400, -200, 200 ) );
85  AddHist( mHists["hcaldY"] = new TH1F( "hcaldY", ";dy (trk - hcal); counts", 400, -200, 200 ) );
86  AddHist( mHists["hcaldYdNFit"] = new TH2F( "hcaldYdNFit", ";dy (trk - hcal); nFit", 400, -200, 200, 10, 0, 10 ) );
87  AddHist( mHists["matchedHcaldY"] = new TH1F( "matchedHcaldY", ";dy (trk - hcal); counts", 400, -200, 200 ) );
88  AddHist( mHists["hcaldR"] = new TH1F( "hcaldR", ";dr (trk - hcal); counts", 400, 0, 400 ) );
89  AddHist( mHists["hcalMindR"] = new TH1F( "hcalMindR", ";dr (trk - hcal); counts", 400, 0, 400 ) );
90  AddHist( mHists["matchedHcaldR"] = new TH1F( "matchedHcaldR", ";dr (trk - hcal); counts", 400, 0, 400 ) );
91  AddHist( mHists["trkEcalX"] = new TH2F( "trkEcalX", ";trkX;ecalX", 300, -150, 150, 300, -150, 150 ) );
92  AddHist( mHists["trkEcalY"] = new TH2F( "trkEcalY", ";trkY;ecalY", 300, -150, 150, 300, -150, 150 ) );
93  AddHist( mHists["trkEcalMinX"] = new TH2F( "trkEcalMinX", ";trkX;ecalX", 300, -150, 150, 300, -150, 150 ) );
94  AddHist( mHists["trkEcalMinY"] = new TH2F( "trkEcalMinY", ";trkY;ecalY", 300, -150, 150, 300, -150, 150 ) );
95  AddHist( mHists["trkHcalX"] = new TH2F( "trkHcalX", ";trkX;hcalX", 300, -150, 150, 300, -150, 150 ) );
96  AddHist( mHists["trkHcalY"] = new TH2F( "trkHcalY", ";trkY;hcalY", 300, -150, 150, 300, -150, 150 ) );
97  AddHist( mHists["trkHcalMinX"] = new TH2F( "trkHcalMinX", ";trkX;hcalX", 300, -150, 150, 300, -150, 150 ) );
98  AddHist( mHists["trkHcalMinY"] = new TH2F( "trkHcalMinY", ";trkY;hcalY", 300, -150, 150, 300, -150, 150 ) );
99 
100  return kStOK;
101 }
102 //________________________________________________________________________
104  LOG_DEBUG << "StFwdAnalysisMaker::Make" << endm;
105  long long itStart = FwdTrackerUtils::nowNanoSecond();
106  if (!mAnalyzeMuDst)
107  ProcessFwdTracks();
108  else
109  ProcessFwdMuTracks();
110  LOG_DEBUG << "Processing Fwd Tracks took: " << (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e6 << " ms" << endm;
111  return kStOK;
112 } // Make
113 //________________________________________________________________________
114 void StFwdAnalysisMaker::Clear(const Option_t *opts) { LOG_DEBUG << "StFwdAnalysisMaker::CLEAR" << endm; }
115 //________________________________________________________________________
116 void StFwdAnalysisMaker::ProcessFwdTracks( ){
117  // This is an example of how to process fwd track collection
118  LOG_DEBUG << "StFwdAnalysisMaker::ProcessFwdTracks" << endm;
119  StEvent *stEvent = static_cast<StEvent *>(GetInputDS("StEvent"));
120  if (!stEvent)
121  return;
122 
123  if (stEvent){
124  StFttCollection *fttCol = stEvent->fttCollection();
125  if (fttCol){
126  LOG_DEBUG << "The Ftt Collection has " << fttCol->numberOfPoints() << " points" << endm;
127  }
128  }
129  StFwdTrackCollection * ftc = stEvent->fwdTrackCollection();
130  if (!ftc) {
131  LOG_DEBUG << "Forward Track Collection is not present" << endm;
132  return;
133  }
134 
135  LOG_DEBUG << "Checking FcsCollection" << endm;
136  StFcsCollection *fcs = stEvent->fcsCollection();
137  if (!fcs) return;
138 
139  StFcsDb *mFcsDb = static_cast<StFcsDb *>(GetDataSet("fcsDb"));
140 
141  size_t fwdMultEcalMatch = 0;
142  size_t fwdMultHcalMatch = 0;
143  size_t fwdMultFST = 0;
144 
145  LOG_INFO << "FwdTrackCollection has: " << ftc->tracks().size() << " tracks" << endm;
146 
147  getHist( "fwdMultAll" )->Fill( ftc->tracks().size() );
148 
149  // Cluster info (independen t of tracks)
150  size_t fwdMultEcalClusters = 0;
151  size_t fwdMultHcalClusters = 0;
152  for ( int iDet = 0; iDet < 4; iDet++ ){
153  for( size_t i = 0; i < fcs->clusters(iDet).size(); i++){
154  StFcsCluster * clu = fcs->clusters(iDet)[i];
155 
156  if ( iDet < 2 ){
157  fwdMultEcalClusters++;
158  getHist( "ecalEnergy" )->Fill( clu->energy() );
159  } else if ( iDet < 4 ){
160  fwdMultHcalClusters++;
161  getHist( "hcalEnergy" )->Fill( clu->energy() );
162  }
163  }
164  }
165 
166  getHist( "fwdMultEcalClusters" )->Fill( fwdMultEcalClusters );
167  getHist( "fwdMultHcalClusters" )->Fill( fwdMultHcalClusters );
168 
169 
170  size_t nGood = 0;
171  size_t nFailed = 0;
172  for ( auto fwdTrack : ftc->tracks() ){
173  if ( !fwdTrack->didFitConvergeFully() ) {
174  nFailed++;
175  continue;
176  }
177  nGood++;
178  LOG_DEBUG << TString::Format("StFwdTrack[ nProjections=%lu, nFTTSeeds=%lu, nFSTSeeds=%lu, mPt=%f ]", fwdTrack->mProjections.size(), fwdTrack->mFTTPoints.size(), fwdTrack->mFSTPoints.size(), fwdTrack->momentum().perp()) << endm;
179  LOG_DEBUG << "track fit momentum " << TString::Format( "(pt=%f, eta=%f, phi=%f)", fwdTrack->momentum().perp(), fwdTrack->momentum().pseudoRapidity(), fwdTrack->momentum().phi() ) << endm;
180  LOG_DEBUG << "StFwdTrack has " << fwdTrack->ecalClusters().size() << " ecal matches" << endm;
181  LOG_DEBUG << "StFwdTrack has " << fwdTrack->hcalClusters().size() << " hcal matches" << endm;
182 
183  getHist("ecalMatchPerTrack")->Fill( fwdTrack->ecalClusters().size() );
184  getHist("hcalMatchPerTrack")->Fill( fwdTrack->hcalClusters().size() );
185 
186  getHist( "nHitsFit" )->Fill( fwdTrack->numberOfFitPoints() );
187 
188  if (fwdTrack->mFSTPoints.size() > 0){
189  fwdMultFST ++;
190  }
191 
192  getHist("eta")->Fill( fwdTrack->momentum().pseudoRapidity() );
193  getHist("phi")->Fill( fwdTrack->momentum().phi() );
194  getHist("pt")->Fill( fwdTrack->momentum().perp() );
195 
196  getHist("charge")->Fill( fwdTrack->charge() );
197 
198  // ecal proj
199  int detId = kFcsWcalId;
200  TVector3 ecalXYZ;
201  TVector3 ecapP;
202 
203  StFwdTrackProjection ecalProj = fwdTrack->getProjectionFor( detId, 0 );
204  StFwdTrackProjection hcalProj = fwdTrack->getProjectionFor( kFcsHcalId, 0 );
205  LOG_DEBUG << "EcalProj z= " << ecalProj.mXYZ.z() << endm;
206  LOG_DEBUG << "HcalProj z= " << hcalProj.mXYZ.z() << endm;
207  LOG_DEBUG << "EcalProj Mom" << TString::Format( "(pt=%f, eta=%f, phi=%f)", ecalProj.mMom.perp(), ecalProj.mMom.pseudoRapidity(), ecalProj.mMom.phi() ) << endm;
208 
209  for ( size_t iEcal = 0; iEcal < fwdTrack->ecalClusters().size(); iEcal++ ){
210  StFcsCluster *clu = fwdTrack->ecalClusters()[iEcal];
211  LOG_DEBUG << "Ecal clu detId = " << clu->detectorId() << endm;
212  getHist("matchedEcalEnergy")->Fill( clu->energy() );
213 
214  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(clu->detectorId(), clu->x(), clu->y());
215  float dx = ecalProj.mXYZ.x() - xyz.x();
216  float dy = ecalProj.mXYZ.y() - xyz.y();
217  float dr = sqrt(dx*dx + dy*dy);
218  getHist("matchedEcaldX")->Fill( dx );
219  getHist("matchedEcaldY")->Fill( dy );
220  getHist("matchedEcaldR")->Fill( dr );
221  }
222 
223  if (ecalProj.mXYZ.z() > 500){
224  double mindR = 999;
225  StFcsCluster * cclu = nullptr; // closet cluster
226  for ( int iDet = 0; iDet < 2; iDet++ ){
227  for( size_t i = 0; i < fcs->clusters(iDet).size(); i++){
228  StFcsCluster * clu = fcs->clusters(iDet)[i];
229 
230  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(clu->detectorId(), clu->x(), clu->y());
231  getHist("ecalXY")->Fill( xyz.x(), xyz.y() );
232 
233  float dx = ecalProj.mXYZ.x() - xyz.x();
234  float dy = ecalProj.mXYZ.y() - xyz.y();
235  float dr = sqrt(dx*dx + dy*dy);
236 
237  if ( fabs(dy) < 25 )
238  getHist( "ecaldX" )->Fill( dx );
239  if ( fabs(dx) < 25 )
240  getHist( "ecaldY" )->Fill( dy );
241  getHist( "ecaldR" )->Fill( dr );
242  if ( dr < mindR ){
243  mindR = dr;
244  cclu = clu;
245  }
246 
247  getHist( "trkEcalX" ) -> Fill( ecalProj.mXYZ.x(), xyz.x() );
248  getHist( "trkEcalY" ) -> Fill( ecalProj.mXYZ.y(), xyz.y() );
249 
250  }
251  }
252  getHist( "ecalMindR" )->Fill( mindR );
253  if (cclu){
254  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(cclu->detectorId(), cclu->x(), cclu->y());
255  getHist( "trkEcalMinX" ) -> Fill( ecalProj.mXYZ.x(), xyz.x() );
256  getHist( "trkEcalMinY" ) -> Fill( ecalProj.mXYZ.y(), xyz.y() );
257  }
258  }
259 
260  if (hcalProj.mXYZ.z() > 500){
261 
262  double mindR = 999;
263  StFcsCluster * cclu = nullptr;
264  for ( int iDet = 2; iDet < 4; iDet++ ){
265  for( size_t i = 0; i < fcs->clusters(iDet).size(); i++){
266  StFcsCluster * clu = fcs->clusters(iDet)[i];
267  if (!clu) continue;
268  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(clu->detectorId(), clu->x(), clu->y());
269  getHist("hcalXY")->Fill( xyz.x(), xyz.y() );
270 
271  float dx = hcalProj.mXYZ.x() - xyz.x();
272  float dy = hcalProj.mXYZ.y() - xyz.y();
273  float dr = sqrt(dx*dx + dy*dy);
274 
275  if ( fabs(dy) < 25 ){
276  getHist( "hcaldX" )->Fill( dx );
277  getHist( "hcaldXdNFit" )->Fill( dx, fwdTrack->numberOfFitPoints() );
278 
279  }
280  if ( fabs(dx) < 25 ){
281  getHist( "hcaldY" )->Fill( dy );
282  getHist( "hcaldYdNFit" )->Fill( dy, fwdTrack->numberOfFitPoints() );
283  }
284  getHist( "hcaldR" )->Fill( dr );
285 
286  if ( dr < mindR ){
287  mindR = dr;
288  cclu = clu;
289  }
290 
291  getHist( "trkHcalX" ) -> Fill( hcalProj.mXYZ.x(), xyz.x() );
292  getHist( "trkHcalY" ) -> Fill( hcalProj.mXYZ.y(), xyz.y() );
293  }
294  }
295  getHist( "hcalMindR" )->Fill( mindR );
296  if (cclu){
297  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(cclu->detectorId(), cclu->x(), cclu->y());
298  getHist( "trkHcalMinX" ) -> Fill( hcalProj.mXYZ.x(), xyz.x() );
299  getHist( "trkHcalMinY" ) -> Fill( hcalProj.mXYZ.y(), xyz.y() );
300  }
301  }
302 
303  if (fwdTrack->ecalClusters().size() > 0)
304  fwdMultEcalMatch++;
305  if (fwdTrack->hcalClusters().size() > 0)
306  fwdMultHcalMatch++;
307 
308  } // Loop ftc->tracks()
309 
310  getHist( "fwdMultGood" )->Fill( nGood );
311  getHist( "fwdMultFailed" )->Fill( nFailed );
312  getHist("fwdMultFST")->Fill( fwdMultFST );
313  getHist("fwdMultHcalMatch")->Fill( fwdMultHcalMatch );
314  getHist("fwdMultEcalMatch")->Fill( fwdMultEcalMatch );
315 
316  LOG_INFO << "Found " << nFailed << " failed track fits out of " << ftc->tracks().size() << endm;
317 } // ProcessFwdTracks
318 
319 //________________________________________________________________________
320 void StFwdAnalysisMaker::ProcessFwdMuTracks( ){
321  // This is an example of how to process fwd track collection
322  LOG_DEBUG << "StFwdAnalysisMaker::ProcessFwdMuTracks" << endm;
323  StMuDstMaker *mMuDstMaker = (StMuDstMaker *)GetMaker("MuDst");
324  if(!mMuDstMaker) {
325  LOG_WARN << " No MuDstMaker ... bye-bye" << endm;
326  return;
327  }
328  StMuDst *mMuDst = mMuDstMaker->muDst();
329  if(!mMuDst) {
330  LOG_WARN << " No MuDst ... bye-bye" << endm;
331  return;
332  }
334  if (!ftc) return;
335 
336  StMuFcsCollection *fcs = mMuDst->muFcsCollection();
337  if (!fcs) return;
338 
339  LOG_INFO << "Number of StMuFwdTracks: " << ftc->numberOfFwdTracks() << endl;
340 
341  StFcsDb *mFcsDb = static_cast<StFcsDb *>(GetDataSet("fcsDb"));
342 
343  size_t fwdMultFST = 0;
344  size_t fwdMultEcalMatch = 0;
345  size_t fwdMultHcalMatch = 0;
346 
347  for ( size_t iTrack = 0; iTrack < ftc->numberOfFwdTracks(); iTrack++ ){
348  StMuFwdTrack * muFwdTrack = ftc->getFwdTrack( iTrack );
349  // LOG_DEBUG << TString::Format("StMuFwdTrack[ nProjections=%lu, nFTTSeeds=%lu, nFSTSeeds=%lu, mPt=%f ]", muFwdTrack->mProjections.size(), muFwdTrack->mFTTPoints.size(), muFwdTrack->mFSTPoints.size(), muFwdTrack->momentum().Pt()) << endm;
350 
351  LOG_DEBUG << "StMuFwdTrack has " << muFwdTrack->mEcalClusters.GetEntries() << " Ecal matched" << endm;
352  LOG_DEBUG << "StMuFwdTrack has " << muFwdTrack->mHcalClusters.GetEntries() << " Hcal matched" << endm;
353 
354  getHist("eta")->Fill( muFwdTrack->momentum().Eta() );
355  getHist("phi")->Fill( muFwdTrack->momentum().Phi() );
356 
357  if (muFwdTrack->mFSTPoints.size() > 0){
358  fwdMultFST ++;
359  }
360 
361  if (muFwdTrack->mEcalClusters.GetEntries() > 0)
362  fwdMultEcalMatch++;
363  if (muFwdTrack->mHcalClusters.GetEntries() > 0)
364  fwdMultHcalMatch++;
365 
366 
367  // ecal proj
368  int detId = kFcsWcalId;
369  TVector3 ecalXYZ;
370  TVector3 ecapP;
371 
372  StMuFwdTrackProjection ecalProj;
373  bool foundEcalProj = muFwdTrack->getProjectionFor( detId, ecalProj, 0 );
374 
375  if (foundEcalProj){
376  for( size_t i = 0; i < fcs->numberOfClusters(); i++){
377  StMuFcsCluster * clu = fcs->getCluster(i);
378 
379  if ( clu->detectorId() > 1 ) continue;
380 
381  if ( clu->energy() < 1 ) continue;
382  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(clu->detectorId(), clu->x(), clu->y());
383 
384  float dx = ecalProj.mXYZ.X() - xyz.x();
385  float dy = ecalProj.mXYZ.Y() - xyz.y();
386  float dr = sqrt(dx*dx + dy*dy);
387 
388  getHist( "ecaldX" )->Fill( dx );
389  getHist( "ecaldY" )->Fill( dy );
390  getHist( "ecaldR" )->Fill( dr );
391 
392  getHist( "trkEcalX" ) -> Fill( ecalProj.mXYZ.X(), xyz.x() );
393 
394  } // i
395  } // foundEcalProj
396 
397 
398  for ( int i = 0; i < muFwdTrack->mEcalClusters.GetEntries(); i++ ){
399  auto c = (StMuFcsCluster*) muFwdTrack->mEcalClusters.At(i);
400  if (!c) continue;
401  getHist("ecalEnergy")->Fill( c->energy() );
402 
403  LOG_DEBUG << "eCal Cluster detId = " << c->detectorId() << endm;
404  StThreeVectorD xyz = mFcsDb->getStarXYZfromColumnRow(c->detectorId(), c->x(), c->y());
405  getHist("ecalXY")->Fill( xyz.x(), xyz.y() );
406 
407  if (foundEcalProj){
408  getHist("matchedEcaldX")->Fill( ecalProj.mXYZ.X() - xyz.x() );
409  }
410  } // i
411 
412  getHist("ecalMatchPerTrack")->Fill( muFwdTrack->mEcalClusters.GetEntries() );
413  getHist("hcalMatchPerTrack")->Fill( muFwdTrack->mHcalClusters.GetEntries() );
414 
415  for ( int i = 0; i < muFwdTrack->mHcalClusters.GetEntries(); i++ ){
416  auto c = (StMuFcsCluster*) muFwdTrack->mHcalClusters.At(i);
417  if (!c) continue;
418  getHist("hcalEnergy")->Fill( c->energy() );
419 
420  getHist("hcalXY")->Fill( c->x(), c->y() );
421  } // i
422  } // iTrack
423 
424  getHist("fwdMult")->Fill( ftc->numberOfFwdTracks() );
425  getHist("fwdMultFST")->Fill( fwdMultFST );
426  getHist("fwdMultHcalMatch")->Fill( fwdMultHcalMatch );
427  getHist("fwdMultEcalMatch")->Fill( fwdMultEcalMatch );
428 }
static StMuFcsCollection * muFcsCollection()
returns pointer to current StMuFcsCollection
Definition: StMuDst.h:395
StMuDst * muDst()
Definition: StMuDstMaker.h:425
std::map< TString, TH1 * > mHists
Map of &lt;name (TString), histogram&gt;
StThreeVectorD getStarXYZfromColumnRow(int det, float col, float row, float FcsZ=-1.0) const
get the STAR frame cooridnates from other way
Definition: StFcsDb.cxx:863
static StMuFwdTrackCollection * muFwdTrackCollection()
returns pointer to current StMuFwdTrackCollection
Definition: StMuDst.h:401
bool mAnalyzeMuDst
Control whether the analysis uses StEvent (default) or MuDst as input.
TH1 * getHist(TString n)
Get the Hist object from the map.
Definition: Stypes.h:40
Definition: Stypes.h:41