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