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"), _etaBins(49), _etaMax(4.9) { 00030 } 00031 00032 00033 /// @name Analysis methods 00034 //@{ 00035 /// Book histograms and initialise projections before the run 00036 void init() { 00037 //All final states. Rapidity range = ATLAS calorimetry. Lowest pT cut = 200 MeV. 00038 const FinalState cnfs2(-_etaMax, _etaMax, 0.2 * GeV); 00039 const FinalState cnfs4(-_etaMax, _etaMax, 0.4 * GeV); 00040 const FinalState cnfs6(-_etaMax, _etaMax, 0.6 * GeV); 00041 const FinalState cnfs8(-_etaMax, _etaMax, 0.8 * GeV); 00042 addProjection(cnfs2, "CNFS2"); 00043 addProjection(cnfs4, "CNFS4"); 00044 addProjection(cnfs6, "CNFS6"); 00045 addProjection(cnfs8, "CNFS8"); 00046 00047 _etaBinSize = (2. * _etaMax)/(double)_etaBins; 00048 00049 //Book histogram 00050 _h_DeltaEtaF_200 = bookHisto1D(1, 1, 1); 00051 _h_DeltaEtaF_400 = bookHisto1D(2, 1, 1); 00052 _h_DeltaEtaF_600 = bookHisto1D(3, 1, 1); 00053 _h_DeltaEtaF_800 = bookHisto1D(4, 1, 1); 00054 } 00055 00056 private: 00057 void fillMap(const FinalState& fs, bool* energyMap, double pTcut) { 00058 // Fill true/false array by iterating over particles and compare their 00059 // pT with pTcut 00060 foreach (const Particle& p, fs.particles(cmpMomByEta)) { 00061 int checkBin = -1; 00062 double checkEta = -_etaMax; 00063 while (1) { 00064 checkEta += _etaBinSize; 00065 ++checkBin; 00066 if (p.eta() < checkEta) { 00067 energyMap[checkBin] = (p.pT() > pTcut * GeV); 00068 break; 00069 } 00070 } 00071 } 00072 } 00073 00074 public: 00075 /// Perform the per-event analysis 00076 void analyze(const Event& event) { 00077 static unsigned int event_count = 0; 00078 ++event_count; 00079 const double weight = event.weight(); 00080 const FinalState& fs2 = applyProjection<FinalState>(event, "CNFS2"); 00081 const FinalState& fs4 = applyProjection<FinalState>(event, "CNFS4"); 00082 const FinalState& fs6 = applyProjection<FinalState>(event, "CNFS6"); 00083 const FinalState& fs8 = applyProjection<FinalState>(event, "CNFS8"); 00084 00085 // Set up Yes/No arrays for energy in each eta bin at each pT cut 00086 bool* energyMap_200 = new bool[_etaBins]; 00087 bool* energyMap_400 = new bool[_etaBins]; 00088 bool* energyMap_600 = new bool[_etaBins]; 00089 bool* energyMap_800 = new bool[_etaBins]; 00090 for (int i = 0; i < _etaBins; ++i) { 00091 energyMap_200[i] = false; 00092 energyMap_400[i] = false; 00093 energyMap_600[i] = false; 00094 energyMap_800[i] = false; 00095 } 00096 00097 // Veto bins based on final state particles > Cut (Where Cut = 200 - 800 MeV pT) 00098 fillMap(fs2, energyMap_200, 0.2); 00099 fillMap(fs4, energyMap_400, 0.4); 00100 fillMap(fs6, energyMap_600, 0.6); 00101 fillMap(fs8, energyMap_800, 0.8); 00102 00103 // Apply gap finding algorithm 00104 // Detector layout follows... 00105 // -Eta [Proton --- DetectorCSide --- DetectorBarrel --- DetectorASide --- Proton] +Eta 00106 bool gapDirectionAt200 = false; //False is gap on C size, True is gap on A side. 00107 double largestEdgeGap_200 = 0.; 00108 double largestEdgeGap_400 = 0.; 00109 double largestEdgeGap_600 = 0.; 00110 double largestEdgeGap_800 = 0.; 00111 00112 for (int E = 200; E <= 800; E += 200) { 00113 double EdgeGapSizeA = -1, EdgeGapSizeC = -1; 00114 bool* energyMap = 0; 00115 switch (E) { 00116 case 200: energyMap = energyMap_200; break; 00117 case 400: energyMap = energyMap_400; break; 00118 case 600: energyMap = energyMap_600; break; 00119 case 800: energyMap = energyMap_800; break; 00120 } 00121 00122 // Look left to right 00123 for (int a = 0; a < _etaBins; ++a) { 00124 if (energyMap[a] == true) { 00125 EdgeGapSizeA = (_etaBinSize * a); 00126 break; 00127 } 00128 } 00129 00130 // And look right to left 00131 for (int c = _etaBins-1; c >= 0; --c) { 00132 if (energyMap[c] == true) { 00133 EdgeGapSizeC = (2 * _etaMax) - (_etaBinSize * (c+1)); 00134 if (fuzzyEquals(EdgeGapSizeC, 4.47035e-08)) EdgeGapSizeC = 0.0; 00135 break; 00136 } 00137 } 00138 // Put your hands on your hips 00139 00140 // Find the largest gap 00141 double largestEdgeGap = 0.; 00142 if (E == 200) { 00143 // If the 200 MeV pass, take the biggest of the two gaps. Make note of which side for higher pT cuts. 00144 largestEdgeGap = std::max(EdgeGapSizeA,EdgeGapSizeC); 00145 gapDirectionAt200 = (EdgeGapSizeA > EdgeGapSizeC); 00146 } else { 00147 // Use the direction from 200 MeV pass, most accurate measure of which side gap is on. 00148 if (gapDirectionAt200) { 00149 largestEdgeGap = EdgeGapSizeA; 00150 } 00151 else largestEdgeGap = EdgeGapSizeC; 00152 } 00153 00154 // Check case of empty detector 00155 if (largestEdgeGap < 0.0) largestEdgeGap = 2.0 * _etaMax; 00156 00157 // Fill bin centre 00158 switch (E) { 00159 case 200: _h_DeltaEtaF_200->fill(largestEdgeGap + _etaBinSize/2., weight); break; 00160 case 400: _h_DeltaEtaF_400->fill(largestEdgeGap + _etaBinSize/2., weight); break; 00161 case 600: _h_DeltaEtaF_600->fill(largestEdgeGap + _etaBinSize/2., weight); break; 00162 case 800: _h_DeltaEtaF_800->fill(largestEdgeGap + _etaBinSize/2., weight); break; 00163 } 00164 00165 if (E == 200) largestEdgeGap_200 = largestEdgeGap; 00166 if (E == 400) largestEdgeGap_400 = largestEdgeGap; 00167 if (E == 600) largestEdgeGap_600 = largestEdgeGap; 00168 if (E == 800) largestEdgeGap_800 = largestEdgeGap; 00169 } 00170 00171 // Algorithm result every 1000 events 00172 if (event_count % 1000 == 0) { 00173 for (int E = 200; E <= 800; E += 200) { 00174 bool* energyMap = 0; 00175 double largestEdgeGap = 0; 00176 switch (E) { 00177 case 200: energyMap = energyMap_200; largestEdgeGap = largestEdgeGap_200; break; 00178 case 400: energyMap = energyMap_400; largestEdgeGap = largestEdgeGap_400; break; 00179 case 600: energyMap = energyMap_600; largestEdgeGap = largestEdgeGap_600; break; 00180 case 800: energyMap = energyMap_800; largestEdgeGap = largestEdgeGap_800; break; 00181 } 00182 MSG_DEBUG("Largest Forward Gap at pT Cut " << E << " MeV=" << largestEdgeGap 00183 << " eta, NFinalState pT > 200 in ATLAS acceptance:" << fs2.particles().size()); 00184 std::string hitPattern = "Detector HitPattern=-4.9["; 00185 for (int a = 0; a < _etaBins; ++a) { 00186 if (energyMap[a] == true) hitPattern += "X"; 00187 else hitPattern += "_"; 00188 } 00189 hitPattern += "]4.9"; 00190 MSG_DEBUG(hitPattern); 00191 std::string gapArrow = " "; 00192 if (!gapDirectionAt200) { 00193 int drawSpaces = (int)(_etaBins - (largestEdgeGap/_etaBinSize) + 0.5); 00194 for (int i = 0; i < drawSpaces; ++i) gapArrow += " "; 00195 } 00196 int drawArrows = (int)((largestEdgeGap/_etaBinSize) + 0.5); 00197 for (int i = 0; i < drawArrows; ++i) gapArrow += "^"; 00198 MSG_DEBUG(gapArrow); 00199 } 00200 } 00201 00202 delete[] energyMap_200; 00203 delete[] energyMap_400; 00204 delete[] energyMap_600; 00205 delete[] energyMap_800; 00206 } 00207 00208 /// Normalise histograms after the run, Scale to cross section 00209 void finalize() { 00210 MSG_DEBUG("Cross Section=" << crossSection() / millibarn << "mb, SumOfWeights=" << sumOfWeights()); 00211 scale(_h_DeltaEtaF_200, (crossSection() / millibarn)/sumOfWeights()); 00212 scale(_h_DeltaEtaF_400, (crossSection() / millibarn)/sumOfWeights()); 00213 scale(_h_DeltaEtaF_600, (crossSection() / millibarn)/sumOfWeights()); 00214 scale(_h_DeltaEtaF_800, (crossSection() / millibarn)/sumOfWeights()); 00215 } 00216 //@} 00217 00218 private: 00219 /// @name Histograms 00220 //@{ 00221 Histo1DPtr _h_DeltaEtaF_200; 00222 Histo1DPtr _h_DeltaEtaF_400; 00223 Histo1DPtr _h_DeltaEtaF_600; 00224 Histo1DPtr _h_DeltaEtaF_800; 00225 //@} 00226 /// @name Private variables 00227 //@{ 00228 int _etaBins; 00229 double _etaMax; 00230 double _etaBinSize; 00231 //@} 00232 }; 00233 00234 // The hook for the plugin system 00235 DECLARE_RIVET_PLUGIN(ATLAS_2012_I1084540); 00236 00237 } Generated on Tue Mar 24 2015 17:35:24 for The Rivet MC analysis system by ![]() |