CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/QCDAnalysis/ChargedHadronSpectra/src/Histograms.cc

Go to the documentation of this file.
00001 #include "QCDAnalysis/ChargedHadronSpectra/interface/Histograms.h"
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 
00004 #include "TROOT.h"
00005 #include "TTree.h"
00006 #include "TFile.h"
00007 #include "TH1F.h"
00008 #include "TH2F.h"
00009 #include "TH3F.h"
00010 
00011 #include <iostream>
00012 
00013 #include <cmath>
00014 using namespace std;
00015 
00016 /*****************************************************************************/
00017 #define nCharges 2
00018 enum Charges
00019 { pos, neg, zero, undefined  };
00020 
00021 const char *chargeName[nCharges + 2] =
00022 { "pos", "neg", "zero", "undefined" };
00023 
00024 /*****************************************************************************/
00025 #define nParticles 21
00026 enum Particles
00027 {
00028   pip, pim, kap, kam,
00029   prp, prm, elp, elm,
00030   hap, ham,
00031   gam, k0s, lam, ala,
00032   rho, kst, aks, phi,
00033   sip, asi,
00034   any
00035 };
00036 
00037 const char *partName[nParticles] =
00038 {
00039   "pip", "pim", "kap", "kam",
00040   "prp", "prm", "elp", "elm",
00041   "hap", "ham",
00042   "gam", "k0s", "lam", "ala",
00043   "rho", "kst", "aks", "phi",
00044   "sip", "asi",
00045   "any"
00046 };
00047 
00048 const int partCharge[nParticles] =
00049 {
00050   pos , neg , pos , neg ,
00051   pos , neg , pos , neg ,
00052   pos , neg ,
00053   zero, zero, zero, zero,
00054   zero, zero, zero, zero,
00055   pos , neg ,
00056   undefined
00057 };
00058 
00059 /*****************************************************************************/
00060 #define nFeedDowns 18
00061 const std::pair<int,int> feedDown[nFeedDowns] =
00062 {
00063   std::pair<int,int>(k0s, pip), std::pair<int,int>(k0s, pim),
00064   std::pair<int,int>(lam, prp), std::pair<int,int>(lam, pim),
00065   std::pair<int,int>(ala, prm), std::pair<int,int>(ala, pip),
00066   std::pair<int,int>(sip, prp), std::pair<int,int>(asi, prm),
00067   std::pair<int,int>(any, pip), std::pair<int,int>(any, pim),
00068   std::pair<int,int>(any, kap), std::pair<int,int>(any, kam),
00069   std::pair<int,int>(any, prp), std::pair<int,int>(any, prm),
00070   std::pair<int,int>(any, elp), std::pair<int,int>(any, elm),
00071   std::pair<int,int>(any, hap), std::pair<int,int>(any, ham)
00072 };
00073 
00074 #define nResonances 4
00075 const std::pair<int,int> resonance[nResonances] =
00076 {
00077  std::pair<int,int>(pip, pim),
00078  std::pair<int,int>(kap, pim),
00079  std::pair<int,int>(pip, kam),
00080  std::pair<int,int>(kap, kam)
00081 };
00082 
00083 /*****************************************************************************/
00084 Histograms::Histograms(const edm::ParameterSet& pset)
00085 {
00086   fillHistograms         = pset.getParameter<bool>("fillHistograms");
00087   fillNtuples            = pset.getParameter<bool>("fillNtuples");
00088 
00089   std::string histoFileLabel = pset.getParameter<std::string>("histoFile");
00090   histoFile = new TFile(histoFileLabel.c_str(),"recreate");
00091 
00092   std::string ntupleFileLabel = pset.getParameter<std::string>("ntupleFile");
00093   ntupleFile = new TFile(ntupleFileLabel.c_str(),"recreate");
00094 
00095 //  resultFile->cd();
00096 }
00097 
00098 /*****************************************************************************/
00099 Histograms::~Histograms()
00100 {
00101 }
00102 
00103 /****************************************************************************/
00104 int Histograms::getParticle(int id)
00105 {
00106   switch(id)
00107   {
00108     case   211 : return pip; break;
00109     case  -211 : return pim; break;
00110 
00111     case   321 : return kap; break;
00112     case  -321 : return kam; break;
00113 
00114     case  2212 : return prp; break;
00115     case -2212 : return prm; break;
00116 
00117     case    11 : return elp; break;
00118     case   -11 : return elm; break;
00119 
00120     //  SimG4Core/Notification/src/G4TrackToParticleID.cc
00121     // "deuteron" = -100;  1p 1n
00122     // "triton"   = -101;  1p 2n
00123     // "alpha"    = -102;  2p 2n
00124     // "He3"      = -104;  2p 1n
00125 
00126     case    22 : return gam; break;
00127     case   310 : return k0s; break;
00128     case  3122 : return lam; break;
00129     case -3122 : return ala; break;
00130 
00131 /*
00132     case   113 : return rho; break;
00133     case   313 : return kst; break;
00134     case  -313 : return aks; break;
00135     case   333 : return phi; break;
00136 */
00137 
00138     case  3222 : return sip; break;
00139     case -3222 : return asi; break;
00140 
00141     default    : return -1;  break;
00142   }
00143 }
00144 
00145 /****************************************************************************/
00146 int Histograms::getCharge(int charge)
00147 { 
00148   if(charge > 0) return pos;
00149             else return neg;
00150 }
00151 
00152 /*****************************************************************************/
00153 void Histograms::declareHistograms()
00154 {
00155   if(fillNtuples)
00156   {
00157     TString leafStr;
00158 
00159     trackTrees.push_back(new TTree("simTrackTree","simTrackTree"));
00160     leafStr = "ids/I:etas/F:pts/F:acc/I:prim/I:nrec/I:ntrkr/I";
00161     trackTrees[0]->Branch("simTrackValues", &simTrackValues, leafStr.Data());
00162 
00163     trackTrees.push_back(new TTree("recTrackTree","recTrackTree"));
00164     leafStr = "charge/I:etar/F:ptr/F:phir/F:zr/F:logpr/F:logde/F:nhitr/I:prim/I:nsim/I:ids/I:parids/I:etas/F:pts/F:ntrkr/I";
00165     trackTrees[1]->Branch("recTrackValues", &recTrackValues, leafStr.Data());
00166 
00167     trackTrees.push_back(new TTree("recVzeroTree","recVzeroTree"));
00168     leafStr = "etar/F:ptr/F:ima/F:rhor/F";
00169     trackTrees[2]->Branch("recVzeroValues", &recVzeroValues, leafStr.Data());
00170 
00171     trackTrees.push_back(new TTree("eventInfoTree","eventInfoTree"));
00172     leafStr = "proc/I:strk/I:ntrkr/I";
00173     trackTrees[3]->Branch("eventInfoValues", &eventInfoValues, leafStr.Data());
00174   }
00175 
00176   if(fillHistograms)
00177   {
00179   // Pt
00180   const double small = 1e-3;
00181   double pt;
00182 
00183   for(pt =  0; pt <  1 - small; pt += 0.05) ptBins.push_back(pt);
00184   for(pt =  1; pt <  2 - small; pt += 0.1 ) ptBins.push_back(pt);
00185   for(pt =  2; pt <  4 - small; pt += 0.2 ) ptBins.push_back(pt);
00186   for(pt =  4; pt <  8 - small; pt += 0.5 ) ptBins.push_back(pt);
00187   for(pt =  8; pt < 16 - small; pt += 1.  ) ptBins.push_back(pt);
00188   for(pt = 16; pt < 32 - small; pt += 2.  ) ptBins.push_back(pt);
00189   for(pt = 32; pt < 64 - small; pt += 4.  ) ptBins.push_back(pt);
00190 
00191   static float ratMin   = 0.5;
00192   static float ratMax   = 1.5;
00193   static float ratWidth = 1./200;
00194 
00195   for(double rat = ratMin; rat < ratMax + ratWidth/2; rat += ratWidth)
00196     ratBins.push_back(rat);
00197 
00199   // Eta (-3,3)
00200   static float etaMin   = -3.0;
00201   static float etaMax   =  3.0;
00202   static float etaWidth =  0.2;
00203 
00204   for(double eta = etaMin; eta < etaMax + etaWidth/2; eta += etaWidth)
00205     etaBins.push_back(eta);
00206 
00207 //  for(double eta = etaMin; eta < etaMax + etaWidth/2; eta += etaWidth/10)
00208   for(double eta = etaMin; eta < etaMax + etaWidth/2; eta += etaWidth/5)
00209     metaBins.push_back(eta);
00210 
00211   static float zMin   = -20.;
00212   static float zMax   =  20.;
00213 //  static float zWidth =  0.1;
00214   static float zWidth =  0.2;
00215 
00216   for(double z = zMin; z < zMax + zWidth/2; z += zWidth)
00217     zBins.push_back(z);
00218 
00220   // Number of recontructed tracks
00221   static float ntrkMin   =  0.5;
00222 // FIXME
00223 //  static float ntrkMax   = 200.;
00224 //  static float ntrkWidth =   5.;
00225   static float ntrkMax   = 1000.;
00226   static float ntrkWidth =   10.;
00227   
00228   for(double ntrk = ntrkMin; ntrk < ntrkMax + ntrkWidth; ntrk += ntrkWidth)
00229     ntrkBins.push_back(ntrk);
00230 
00231 
00232   char histName[256];
00233 
00235   // EventInfo
00236   sprintf(histName,"heve");
00237   heve.push_back(new TH1F(histName,histName, 200, -0.5,199.5));
00238 
00239   sprintf(histName,"hsdx");
00240   heve.push_back(new TH1F(histName,histName, 200, -0.5,199.5));
00241 
00242   sprintf(histName,"hddx");
00243   heve.push_back(new TH1F(histName,histName, 200, -0.5,199.5));
00244 
00245   sprintf(histName,"hndx");
00246   heve.push_back(new TH1F(histName,histName, 200, -0.5,199.5));
00247 
00248   sprintf(histName,"hder");
00249   hder.push_back(new TH2F(histName,histName, 200, -0.5,199.5,
00250                                              200, -0.5,199.5));
00251 
00253   // SimTrack
00254   for(int part = pip; part <= ala; part++)
00255   {
00256     // simulated
00257     sprintf(histName,"hsim_%s", partName[part]);
00258     hsim.push_back(new TH3F(histName,histName,
00259                             etaBins.size()-1,  &etaBins[0],
00260                              ptBins.size()-1,   &ptBins[0],
00261                            ntrkBins.size()-1, &ntrkBins[0]));
00262    
00263     // accepted 
00264     sprintf(histName,"hacc_%s", partName[part]);
00265     hacc.push_back(new TH3F(histName,histName,
00266                             etaBins.size()-1,  &etaBins[0],
00267                              ptBins.size()-1,   &ptBins[0],
00268                            ntrkBins.size()-1, &ntrkBins[0]));
00269   
00270     // reconstructed/efficiency
00271     sprintf(histName,"href_%s", partName[part]);
00272     href.push_back(new TH3F(histName,histName,
00273                             etaBins.size()-1,  &etaBins[0],
00274                              ptBins.size()-1,   &ptBins[0],
00275                            ntrkBins.size()-1, &ntrkBins[0]));
00276   
00277     // multiply reconstructed
00278     sprintf(histName,"hmul_%s", partName[part]);
00279     hmul.push_back(new TH3F(histName,histName,
00280                             etaBins.size()-1,  &etaBins[0],
00281                              ptBins.size()-1,   &ptBins[0],
00282                            ntrkBins.size()-1, &ntrkBins[0]));
00283   }
00284 
00286   // RecTrack
00287   for(int charge = 0; charge < nCharges; charge++)
00288   {
00289     sprintf(histName,"hall_%s",chargeName[charge]);
00290     hall.push_back(new TH3F(histName,histName,
00291                             etaBins.size()-1,  &etaBins[0],
00292                              ptBins.size()-1,   &ptBins[0],
00293                            ntrkBins.size()-1, &ntrkBins[0]));
00294 
00295     sprintf(histName,"hdac_%s",chargeName[charge]);
00296     hdac.push_back(new TH2F(histName,histName, 
00297                            metaBins.size()-1, &metaBins[0],
00298                               zBins.size()-1,    &zBins[0]));
00299 
00301     // RecTrack -- FakeRate
00302     sprintf(histName,"hfak_%s",chargeName[charge]);
00303     hfak.push_back(new TH3F(histName,histName,
00304                             etaBins.size()-1,  &etaBins[0],
00305                              ptBins.size()-1,   &ptBins[0],
00306                            ntrkBins.size()-1, &ntrkBins[0]));
00307   }
00308 
00310   // RecTrack -- Resolution, bias
00311   for(int part = pip; part <= ala; part++)
00312   {
00313     // value
00314     sprintf(histName,"hvpt_%s",partName[part]);
00315     hvpt.push_back(new TH3F(histName,histName,
00316                             etaBins.size()-1, &etaBins[0],
00317                              ptBins.size()-1,  &ptBins[0],
00318                              ptBins.size()-1,  &ptBins[0]));
00319 
00320     // ratio
00321     sprintf(histName,"hrpt_%s",partName[part]);
00322     hrpt.push_back(new TH3F(histName,histName,
00323                             etaBins.size()-1, &etaBins[0],
00324                              ptBins.size()-1,  &ptBins[0],
00325                             ratBins.size()-1, &ratBins[0]));
00326 
00327     sprintf(histName,"hsp0_%s",partName[part]);
00328     hsp0.push_back(new TH2F(histName,histName,
00329                             etaBins.size()-1, &etaBins[0],
00330                              ptBins.size()-1,  &ptBins[0]));
00331 
00332     sprintf(histName,"hsp1_%s",partName[part]);
00333     hsp1.push_back(new TH2F(histName,histName,
00334                             etaBins.size()-1, &etaBins[0],
00335                              ptBins.size()-1,  &ptBins[0]));
00336 
00337     sprintf(histName,"hsp2_%s",partName[part]);
00338     hsp2.push_back(new TH2F(histName,histName,
00339                             etaBins.size()-1, &etaBins[0],
00340                              ptBins.size()-1,  &ptBins[0]));
00341   }
00342 
00344   // RecTrack -- FeedDown
00345   for(int k = 0; k < nFeedDowns; k++)
00346   {
00347     sprintf(histName,"hpro_%s_%s", partName[feedDown[k].first], // produced
00348                                    partName[feedDown[k].second]);
00349     hpro.push_back(new TH2F(histName,histName,
00350                             etaBins.size()-1, &etaBins[0],
00351                              ptBins.size()-1,  &ptBins[0]));
00352 
00353     sprintf(histName,"hdec_%s_%s", partName[feedDown[k].first], // decay
00354                                    partName[feedDown[k].second]);
00355     hdec.push_back(new TH2F(histName,histName,
00356                             etaBins.size()-1, &etaBins[0],
00357                              ptBins.size()-1,  &ptBins[0]));
00358   }
00359 
00361   // EnergyLoss
00362   static float lpMin   = -3; // 50 MeV
00363   static float lpMax   =  2; // 7.4 GeV
00364   static float lpWidth = (lpMax - lpMin)/100;
00365   for(double lp = lpMin; lp < lpMax + lpWidth/2; lp += lpWidth)
00366     lpBins.push_back(lp);
00367 
00368   static float ldeMin   = log(1);
00369   static float ldeMax   = log(100);
00370   static float ldeWidth = (ldeMax - ldeMin)/250;
00371   for(double lde = ldeMin; lde < ldeMax + ldeWidth/2; lde += ldeWidth)
00372     ldeBins.push_back(lde);
00373 
00374   for(double nhit = -0.5; nhit < 50; nhit += 1)
00375     nhitBins.push_back(nhit);
00376 
00377   for(int charge = 0; charge < nCharges; charge++)
00378   {
00379     // All hits
00380     // dE/dx
00381     sprintf(histName,"helo_%s", chargeName[charge]);
00382     helo.push_back(new TH3F(histName,histName,
00383                            etaBins.size()-1, &etaBins[0],
00384                             ptBins.size()-1,  &ptBins[0],
00385                            ldeBins.size()-1, &ldeBins[0]));
00386 
00387     // Number of hits used
00388     sprintf(histName,"hnhi_%s", chargeName[charge]);
00389     hnhi.push_back(new TH3F(histName,histName,
00390                            etaBins.size()-1,  &etaBins[0],
00391                             ptBins.size()-1,   &ptBins[0],
00392                           nhitBins.size()-1, &nhitBins[0]));
00393 
00394     // Demo plot
00395     sprintf(histName,"held_%s", chargeName[charge]);
00396     held.push_back(new TH2F(histName,histName,
00397                             lpBins.size()-1,  &lpBins[0],
00398                            ldeBins.size()-1, &ldeBins[0])); 
00399   }
00400 
00401 /*
00402   for(int charge = 0; charge < nCharges; charge++)
00403   {
00404     // Strip hits
00405     // dE/dx
00406     sprintf(histName,"selo_%s", chargeName[charge]);
00407     selo.push_back(new TH3F(histName,histName,
00408                            etaBins.size()-1, &etaBins[0],
00409                             ptBins.size()-1,  &ptBins[0],
00410                            ldeBins.size()-1, &ldeBins[0]));
00411 
00412     // Number of hits used
00413     sprintf(histName,"snhi_%s", chargeName[charge]);
00414     snhi.push_back(new TH3F(histName,histName,
00415                            etaBins.size()-1,  &etaBins[0],
00416                             ptBins.size()-1,   &ptBins[0],
00417                           nhitBins.size()-1, &nhitBins[0]));
00418     
00419     // Demo plot
00420     sprintf(histName,"seld_%s", chargeName[charge]);
00421     seld.push_back(new TH2F(histName,histName,
00422                             lpBins.size()-1,  &lpBins[0],
00423                            ldeBins.size()-1, &ldeBins[0]));
00424   }
00425 */
00426 
00428   // Invariant mass
00429   static float rhoMin   = 0.;
00430   static float rhoMax   = 5.;
00431   static float rhoWidth = 0.2; 
00432   for(double rho_ = rhoMin; rho_ < rhoMax + rhoWidth/2; rho_ += rhoWidth)
00433     rhoBins.push_back(rho_);
00434 
00435 
00436   for(int part = gam; part <= phi; part++)
00437   {
00438     float imMin = 0;
00439     float imMax = 0;
00440     float imWidth = 0;
00441 
00442     if(part == gam) { imMin = 0.0; imMax = 0.2; imWidth = 0.005; }
00443     if(part == k0s) { imMin = 0.3; imMax = 0.7; imWidth = 0.005; }
00444     if(part == lam ||
00445        part == ala) { imMin = 1.0; imMax = 1.3; imWidth = 0.002; }
00446 
00447     if(part == rho) { imMin = 0.2; imMax = 1.2; imWidth = 0.010; }
00448     if(part == kst ||
00449        part == aks) { imMin = 0.6; imMax = 1.6; imWidth = 0.010; }
00450     if(part == phi) { imMin = 0.9; imMax = 1.1; imWidth = 0.002; }
00451 
00452     std::vector<double> imBins;
00453     double im;
00454     for(im = imMin; im < imMax + imWidth/2; im += imWidth)
00455       imBins.push_back(im);
00456 
00457     if(imWidth > 0)
00458     {
00459       sprintf(histName,"hima_%s", partName[part]);
00460       hima.push_back(new TH3F(histName,histName,
00461                               etaBins.size()-1, &etaBins[0],
00462                                ptBins.size()-1,  &ptBins[0],
00463                                imBins.size()-1,  &imBins[0]));
00464 
00465       if(part >= rho && part <= phi)
00466       {
00467       sprintf(histName,"himp_%s", partName[part]);
00468       hima.push_back(new TH3F(histName,histName,
00469                               etaBins.size()-1, &etaBins[0],
00470                                ptBins.size()-1,  &ptBins[0],
00471                                imBins.size()-1,  &imBins[0]));
00472 
00473       sprintf(histName,"himm_%s", partName[part]);
00474       hima.push_back(new TH3F(histName,histName,
00475                               etaBins.size()-1, &etaBins[0],
00476                                ptBins.size()-1,  &ptBins[0],
00477                                imBins.size()-1,  &imBins[0]));
00478 
00479       sprintf(histName,"himx_%s", partName[part]);
00480       hima.push_back(new TH3F(histName,histName,
00481                               etaBins.size()-1, &etaBins[0],
00482                                ptBins.size()-1,  &ptBins[0],
00483                                imBins.size()-1,  &imBins[0]));
00484       }
00485 
00486       sprintf(histName,"hrho_%s", partName[part]);
00487       hrho.push_back(new TH3F(histName,histName,
00488                               rhoBins.size()-1, &rhoBins[0],
00489                                ptBins.size()-1,  &ptBins[0],
00490                                imBins.size()-1,  &imBins[0]));
00491     }
00492   }
00493   }
00494 }
00495 
00496 /****************************************************************************/
00497 void Histograms::fillEventInfo(int proc, int strk, int ntrk)
00498 {
00499   if(fillNtuples)
00500   {
00501     EventInfo_t e;
00502     e.proc  = proc;
00503     e.strk  = strk;
00504     e.ntrkr = ntrk;
00505 
00506     eventInfoValues = e;
00507 
00508     trackTrees[3]->Fill();
00509   }
00510   
00511   if(fillHistograms)
00512   {
00513   heve[0]->Fill(ntrk);
00514 
00515   if(proc == 92 || proc == 93)
00516     heve[1]->Fill(ntrk); // hsdx
00517 
00518   if(proc == 94)
00519     heve[2]->Fill(ntrk); // hddx
00520 
00521   if(!(proc == 92 || proc == 93 || proc == 94))
00522     heve[3]->Fill(ntrk); // hndx
00523 
00524   // For multiplicity, detector response matrix
00525   hder[0]->Fill(strk,ntrk);
00526   }
00527 }
00528 
00529 /****************************************************************************/
00530 void Histograms::fillSimHistograms(const SimTrack_t & s)
00531 {
00532   if(fillNtuples)
00533   { 
00534     if(s.prim)
00535     {
00536       simTrackValues = s;
00537       trackTrees[0]->Fill();
00538     }
00539   }
00540 
00541   if(fillHistograms)
00542   {
00543     int part = getParticle(s.ids);   
00544 
00545     if(pip <= part && part <= ala && s.prim)
00546     {
00547                      hsim[part]->Fill(s.etas, s.pts, s.ntrkr);
00548       if(s.acc)      hacc[part]->Fill(s.etas, s.pts, s.ntrkr);
00549 if(s.acc)
00550 {
00551       if(s.nrec > 0) href[part]->Fill(s.etas, s.pts, s.ntrkr);
00552       if(s.nrec > 1) hmul[part]->Fill(s.etas, s.pts, s.ntrkr);
00553 }
00554   
00555       if(partCharge[part] == pos || partCharge[part] == neg)
00556       {
00557         if(partCharge[part] == pos) part = hap;
00558         if(partCharge[part] == neg) part = ham;
00559   
00560                        hsim[part]->Fill(s.etas, s.pts, s.ntrkr);
00561         if(s.acc)      hacc[part]->Fill(s.etas, s.pts, s.ntrkr);
00562 if(s.acc)
00563 {
00564         if(s.nrec > 0) href[part]->Fill(s.etas, s.pts, s.ntrkr);
00565         if(s.nrec > 1) hmul[part]->Fill(s.etas, s.pts, s.ntrkr);
00566 }
00567       }
00568     }
00569   }
00570 }
00571 
00572 /****************************************************************************/
00573 void Histograms::fillRecHistograms(const RecTrack_t & r)
00574 {
00575   if(fillNtuples)
00576   {
00577     if(r.prim)
00578     {
00579       recTrackValues = r;
00580       trackTrees[1]->Fill();
00581     }
00582   }
00583 
00584   if(fillHistograms)
00585   {
00586   int charge = getCharge(r.charge);
00587   double p = exp(r.logpr);
00588 
00589   if(r.prim)
00590   {
00591     hall[charge]->Fill(r.etar, r.ptr, r.ntrkr);
00592 
00593     if(r.ptr > 0.3) // !!!
00594     hdac[charge]->Fill(r.etar, r.zr);
00595 
00596     if(r.nsim == 0)
00597     hfak[charge]->Fill(r.etar, r.ptr, r.ntrkr);
00598     
00599     if(r.nsim == 1)
00600     {
00601       int part = getParticle(r.ids);
00602       int moth = getParticle(r.parids);
00603 
00604       if(pip <= part && part <= ala)
00605       {
00606         hvpt[part]->Fill(r.etas, r.pts, r.ptr      ); // value
00607         hrpt[part]->Fill(r.etas, r.pts, r.ptr/r.pts); // ratio 
00608 
00609         hsp0[part]->Fill(r.etas, r.pts);      // sum p^0
00610         hsp1[part]->Fill(r.etas, r.pts, p);   // sum p^1
00611         hsp2[part]->Fill(r.etas, r.pts, p*p); // sum p^2
00612       }
00613 
00614       for(int k = 0; k < nFeedDowns; k++)
00615         if(part == feedDown[k].second) // daughter same
00616         {
00617           hpro[k]->Fill(r.etar, r.ptr);
00618 
00619           if((r.parids != 0 && feedDown[k].first == any) ||
00620                                feedDown[k].first == moth) 
00621           hdec[k]->Fill(r.etar, r.ptr);
00622         }
00623 
00624       if(partCharge[part] == pos || partCharge[part] == neg)
00625       {
00626         if(partCharge[part] == pos) part = hap;
00627         if(partCharge[part] == neg) part = ham;
00628 
00629         hvpt[part]->Fill(r.etas, r.pts, r.ptr      ); // value
00630         hrpt[part]->Fill(r.etas, r.pts, r.ptr/r.pts); // ratio
00631 
00632         hsp0[part]->Fill(r.etas, r.pts);      // sum p^0
00633         hsp1[part]->Fill(r.etas, r.pts, p);   // sum p^1
00634         hsp2[part]->Fill(r.etas, r.pts, p*p); // sum p^2
00635 
00636         for(int k = 0; k < nFeedDowns; k++)
00637         if(part == feedDown[k].second) // daughter same
00638         {
00639           hpro[k]->Fill(r.etar, r.ptr);
00640 
00641           if((r.parids != 0 && feedDown[k].first == any) ||
00642                                feedDown[k].first == moth)
00643           hdec[k]->Fill(r.etar, r.ptr);
00644         }
00645       }
00646     }
00647 
00648     // All hits
00649     helo[charge]->Fill(r.etar, r.ptr, r.logde);
00650     hnhi[charge]->Fill(r.etar, r.ptr, r.nhitr);
00651     held[charge]->Fill(r.logpr, r.logde);
00652 
00653     // Strip hits
00654 /*
00655     selo[charge]->Fill(r.etar, r.ptr, r.logde_strip);
00656     snhi[charge]->Fill(r.etar, r.ptr, r.nhitr_strip);
00657     seld[charge]->Fill(r.logpr, r.logde_strip);
00658 */
00659   }
00660   }
00661 }
00662 
00663 /****************************************************************************/
00664 void Histograms::fillVzeroHistograms(const RecVzero_t & v, int part)
00665 {
00666   if(fillNtuples)
00667   {
00668     recVzeroValues = v;
00669     trackTrees[2]->Fill();
00670   }
00671 
00672   if(fillHistograms)
00673     hima[part]->Fill(v.etar, v.ptr, v.ima);
00674 }
00675 
00676 /****************************************************************************/
00677 void Histograms::writeHistograms()
00678 {
00679   typedef std::vector<TH3F *>::const_iterator H3I;
00680   typedef std::vector<TH2F *>::const_iterator H2I;
00681   typedef std::vector<TH1F *>::const_iterator H1I;
00682   typedef std::vector<TTree *>::const_iterator TI;
00683 
00684   if(fillHistograms)
00685   {
00686   histoFile->cd();
00687   for(H1I h = heve.begin(); h!= heve.end(); h++) (*h)->Write();
00688   for(H2I h = hder.begin(); h!= hder.end(); h++) (*h)->Write();
00689 
00690   for(H3I h = hsim.begin(); h!= hsim.end(); h++) (*h)->Write();
00691   for(H3I h = hacc.begin(); h!= hacc.end(); h++) (*h)->Write();
00692   for(H3I h = href.begin(); h!= href.end(); h++) (*h)->Write();
00693   for(H3I h = hmul.begin(); h!= hmul.end(); h++) (*h)->Write();
00694 
00695   for(H3I h = hall.begin(); h!= hall.end(); h++) (*h)->Write();
00696   for(H2I h = hdac.begin(); h!= hdac.end(); h++) (*h)->Write();
00697 
00698   for(H3I h = hvpt.begin(); h!= hvpt.end(); h++) (*h)->Write();
00699   for(H3I h = hrpt.begin(); h!= hrpt.end(); h++) (*h)->Write();
00700 
00701   for(H2I h = hsp0.begin(); h!= hsp0.end(); h++) (*h)->Write();
00702   for(H2I h = hsp1.begin(); h!= hsp1.end(); h++) (*h)->Write();
00703   for(H2I h = hsp2.begin(); h!= hsp2.end(); h++) (*h)->Write();
00704 
00705   for(H3I h = hfak.begin(); h!= hfak.end(); h++) (*h)->Write();
00706 
00707   for(H2I h = hpro.begin(); h!= hpro.end(); h++) (*h)->Write();
00708   for(H2I h = hdec.begin(); h!= hdec.end(); h++) (*h)->Write();
00709 
00710   for(H3I h = helo.begin(); h!= helo.end(); h++) (*h)->Write();
00711   for(H3I h = hnhi.begin(); h!= hnhi.end(); h++) (*h)->Write();
00712   for(H2I h = held.begin(); h!= held.end(); h++) (*h)->Write();
00713 
00714   for(H3I h = hima.begin(); h!= hima.end(); h++) (*h)->Write();
00715   for(H3I h = hrho.begin(); h!= hrho.end(); h++) (*h)->Write();
00716   histoFile->Close();
00717   }
00718 
00719   if(fillNtuples)
00720   {
00721   ntupleFile->cd();
00722   for(TI t = trackTrees.begin(); t!= trackTrees.end(); t++) (*t)->Write();
00723   ntupleFile->Close();
00724   }
00725 }
00726