rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1084540.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 /**
00003  * @name ATLAS Diffractive Gaps Rivet Analysis
00004  * @author Tim Martin, tim.martin@cern.ch
00005  * @version 1.0
00006  * @date 16/01/2012
00007  * @see http://arxiv.org/abs/1201.2808
00008  * @note pp, sqrt(s) = 7 TeV
00009  * @note Rapidity gap finding algorithm designed to compliment
00010  * the ATLAS detector acceptance. Forward rapidity gap sizes
00011  * are calculated for each event, considering all stable
00012  * particles above pT cut values 200, 400, 600 and 800 MeV in
00013  * turn. A forward rapidity gap is defined to be the largest
00014  * continuous region stretching inward from either edge of the
00015  * detector at eta = |4.9| which contains zero particles above
00016  * pT Cut. Soft diffractive topologies are isolated at large
00017  * gap sizes.
00018  *
00019  */
00020 #include "Rivet/Analysis.hh"
00021 #include "Rivet/Projections/FinalState.hh"
00022 
00023 namespace Rivet {
00024 
00025 
00026   class ATLAS_2012_I1084540 : public Analysis {
00027   public:
00028 
00029     ATLAS_2012_I1084540() : Analysis("ATLAS_2012_I1084540") {}
00030 
00031 
00032     /// @name Analysis methods
00033     //@{
00034     /// Book histograms and initialise projections before the run
00035     void init() {
00036       //All final states. Rapidity range = ATLAS calorimetry. Lowest pT cut = 200 MeV.
00037       const FinalState cnfs2(-_etaMax, _etaMax, 0.2 * GeV);
00038       const FinalState cnfs4(-_etaMax, _etaMax, 0.4 * GeV);
00039       const FinalState cnfs6(-_etaMax, _etaMax, 0.6 * GeV);
00040       const FinalState cnfs8(-_etaMax, _etaMax, 0.8 * GeV);
00041       declare(cnfs2, "CNFS2");
00042       declare(cnfs4, "CNFS4");
00043       declare(cnfs6, "CNFS6");
00044       declare(cnfs8, "CNFS8");
00045 
00046       _etaBinSize = (2. * _etaMax)/(double)_etaBins;
00047 
00048       //Book histogram
00049       _h_DeltaEtaF_200 = bookHisto1D(1, 1, 1);
00050       _h_DeltaEtaF_400 = bookHisto1D(2, 1, 1);
00051       _h_DeltaEtaF_600 = bookHisto1D(3, 1, 1);
00052       _h_DeltaEtaF_800 = bookHisto1D(4, 1, 1);
00053     }
00054 
00055   private:
00056     void fillMap(const FinalState& fs, bool* energyMap, double pTcut) {
00057       // Fill true/false array by iterating over particles and compare their
00058       // pT with pTcut
00059       foreach (const Particle& p, fs.particles(cmpMomByEta)) {
00060         int checkBin = -1;
00061         double checkEta = -_etaMax;
00062         while (1) {
00063           checkEta += _etaBinSize;
00064           ++checkBin;
00065           if (p.eta() < checkEta) {
00066             energyMap[checkBin] = (p.pT() > pTcut * GeV);
00067             break;
00068           }
00069         }
00070       }
00071     }
00072 
00073   public:
00074     /// Perform the per-event analysis
00075     void analyze(const Event& event) {
00076       static unsigned int event_count = 0;
00077       ++event_count;
00078       const double weight = event.weight();
00079       const FinalState& fs2 = apply<FinalState>(event, "CNFS2");
00080       const FinalState& fs4 = apply<FinalState>(event, "CNFS4");
00081       const FinalState& fs6 = apply<FinalState>(event, "CNFS6");
00082       const FinalState& fs8 = apply<FinalState>(event, "CNFS8");
00083 
00084       // Set up Yes/No arrays for energy in each eta bin at each pT cut
00085       bool energyMap_200[_etaBins];
00086       bool energyMap_400[_etaBins];
00087       bool energyMap_600[_etaBins];
00088       bool energyMap_800[_etaBins];
00089       for (int i = 0; i < _etaBins; ++i) {
00090         energyMap_200[i] = false;
00091         energyMap_400[i] = false;
00092         energyMap_600[i] = false;
00093         energyMap_800[i] = false;
00094       }
00095 
00096       // Veto bins based on final state particles > Cut (Where Cut = 200 - 800 MeV pT)
00097       fillMap(fs2, energyMap_200, 0.2);
00098       fillMap(fs4, energyMap_400, 0.4);
00099       fillMap(fs6, energyMap_600, 0.6);
00100       fillMap(fs8, energyMap_800, 0.8);
00101 
00102       // Apply gap finding algorithm
00103       // Detector layout follows...
00104       // -Eta [Proton  ---  DetectorCSide  ---  DetectorBarrel  ---  DetectorASide  ---  Proton] +Eta
00105       bool gapDirectionAt200 = false; //False is gap on C size, True is gap on A side.
00106       double largestEdgeGap_200 = 0.;
00107       double largestEdgeGap_400 = 0.;
00108       double largestEdgeGap_600 = 0.;
00109       double largestEdgeGap_800 = 0.;
00110 
00111       for (int E = 200; E <= 800; E += 200) {
00112         double EdgeGapSizeA = -1, EdgeGapSizeC = -1;
00113         bool* energyMap = 0;
00114         switch (E) {
00115           case 200: energyMap = energyMap_200; break;
00116           case 400: energyMap = energyMap_400; break;
00117           case 600: energyMap = energyMap_600; break;
00118           case 800: energyMap = energyMap_800; break;
00119         }
00120 
00121         // Look left to right
00122         for (int a = 0; a < _etaBins; ++a) {
00123           if (energyMap[a] == true) {
00124             EdgeGapSizeA = (_etaBinSize * a);
00125             break;
00126           }
00127         }
00128 
00129         // And look right to left
00130         for (int c = _etaBins-1; c >= 0; --c) {
00131           if (energyMap[c] == true) {
00132             EdgeGapSizeC = (2 * _etaMax) - (_etaBinSize * (c+1));
00133             if (fuzzyEquals(EdgeGapSizeC, 4.47035e-08)) EdgeGapSizeC = 0.0;
00134             break;
00135           }
00136         }
00137         // Put your hands on your hips
00138 
00139         // Find the largest gap
00140         double largestEdgeGap = 0.;
00141         if (E == 200) {
00142           // If the 200 MeV pass, take the biggest of the two gaps. Make note of which side for higher pT cuts.
00143           largestEdgeGap = std::max(EdgeGapSizeA,EdgeGapSizeC);
00144           gapDirectionAt200 = (EdgeGapSizeA > EdgeGapSizeC);
00145         } else {
00146           // Use the direction from 200 MeV pass, most accurate measure of which side gap is on.
00147           if (gapDirectionAt200) {
00148             largestEdgeGap = EdgeGapSizeA;
00149           }
00150           else largestEdgeGap = EdgeGapSizeC;
00151         }
00152 
00153         // Check case of empty detector
00154         if (largestEdgeGap < 0.0) largestEdgeGap = 2.0 * _etaMax;
00155 
00156         // Fill bin centre
00157         switch (E) {
00158           case 200: _h_DeltaEtaF_200->fill(largestEdgeGap + _etaBinSize/2., weight); break;
00159           case 400: _h_DeltaEtaF_400->fill(largestEdgeGap + _etaBinSize/2., weight); break;
00160           case 600: _h_DeltaEtaF_600->fill(largestEdgeGap + _etaBinSize/2., weight); break;
00161           case 800: _h_DeltaEtaF_800->fill(largestEdgeGap + _etaBinSize/2., weight); break;
00162         }
00163 
00164         if (E == 200) largestEdgeGap_200 = largestEdgeGap;
00165         if (E == 400) largestEdgeGap_400 = largestEdgeGap;
00166         if (E == 600) largestEdgeGap_600 = largestEdgeGap;
00167         if (E == 800) largestEdgeGap_800 = largestEdgeGap;
00168       }
00169 
00170       // Algorithm result every 1000 events
00171       if (event_count % 1000 == 0) {
00172         for (int E = 200; E <= 800; E += 200) {
00173           bool* energyMap = 0;
00174           double largestEdgeGap = 0;
00175           switch (E) {
00176             case 200: energyMap = energyMap_200; largestEdgeGap = largestEdgeGap_200; break;
00177             case 400: energyMap = energyMap_400; largestEdgeGap = largestEdgeGap_400; break;
00178             case 600: energyMap = energyMap_600; largestEdgeGap = largestEdgeGap_600; break;
00179             case 800: energyMap = energyMap_800; largestEdgeGap = largestEdgeGap_800; break;
00180           }
00181           MSG_DEBUG("Largest Forward Gap at pT Cut " << E << " MeV=" << largestEdgeGap
00182             << " eta, NFinalState pT > 200 in ATLAS acceptance:" << fs2.particles().size());
00183           std::string hitPattern = "Detector HitPattern=-4.9[";
00184           for (int a = 0; a < _etaBins; ++a) {
00185             if (energyMap[a] == true) hitPattern += "X";
00186             else hitPattern += "_";
00187           }
00188           hitPattern += "]4.9";
00189           MSG_DEBUG(hitPattern);
00190           std::string gapArrow = "                         ";
00191             if (!gapDirectionAt200) {
00192             int drawSpaces = (int)(_etaBins - (largestEdgeGap/_etaBinSize) + 0.5);
00193             for (int i = 0; i < drawSpaces; ++i) gapArrow += " ";
00194           }
00195           int drawArrows = (int)((largestEdgeGap/_etaBinSize) + 0.5);
00196           for (int i = 0; i < drawArrows; ++i) gapArrow += "^";
00197           MSG_DEBUG(gapArrow);
00198         }
00199       }
00200     }
00201 
00202     /// Normalise histograms after the run, Scale to cross section
00203     void finalize() {
00204       MSG_DEBUG("Cross Section=" << crossSection() / millibarn << "mb, SumOfWeights=" << sumOfWeights());
00205       scale(_h_DeltaEtaF_200, (crossSection() / millibarn)/sumOfWeights());
00206       scale(_h_DeltaEtaF_400, (crossSection() / millibarn)/sumOfWeights());
00207       scale(_h_DeltaEtaF_600, (crossSection() / millibarn)/sumOfWeights());
00208       scale(_h_DeltaEtaF_800, (crossSection() / millibarn)/sumOfWeights());
00209     }
00210     //@}
00211 
00212   private:
00213     /// @name Histograms
00214     //@{
00215     Histo1DPtr _h_DeltaEtaF_200;
00216     Histo1DPtr _h_DeltaEtaF_400;
00217     Histo1DPtr _h_DeltaEtaF_600;
00218     Histo1DPtr _h_DeltaEtaF_800;
00219     //@}
00220     /// @name Private variables
00221     //@{
00222     static constexpr int _etaBins = 49;
00223     static constexpr double _etaMax = 4.9;
00224     double _etaBinSize;
00225     //@}
00226   };
00227 
00228   // The hook for the plugin system
00229   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1084540);
00230 
00231 }