Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "FWCore/Framework/interface/MakerMacros.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012
00013 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00014 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00015 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
00016
00017 #include "SimGeneral/PileupInformation/interface/PileupInformation.h"
00018 #include "SimDataFormats/PileupSummaryInfo/interface/PileupMixingContent.h"
00019
00020 PileupInformation::PileupInformation(const edm::ParameterSet & config)
00021 {
00022
00023
00024 pTcut_1_ = 0.1;
00025 pTcut_2_ = 0.5;
00026 distanceCut_ = config.getParameter<double>("vertexDistanceCut");
00027 volumeRadius_ = config.getParameter<double>("volumeRadius");
00028 volumeZ_ = config.getParameter<double>("volumeZ");
00029 pTcut_1_ = config.getParameter<double>("pTcut_1");
00030 pTcut_2_ = config.getParameter<double>("pTcut_2");
00031
00032 PileupInfoLabel_ = config.getParameter<std::string>("PileupMixingLabel");
00033
00034 trackingTruth_ = config.getParameter<std::string>("TrackingParticlesLabel");
00035
00036 simHitLabel_ = config.getParameter<std::string>("simHitLabel");
00037
00038 MessageCategory_ = "PileupInformation";
00039
00040 edm::LogInfo (MessageCategory_) << "Setting up PileupInformation";
00041 edm::LogInfo (MessageCategory_) << "Vertex distance cut set to " << distanceCut_ << " mm";
00042 edm::LogInfo (MessageCategory_) << "Volume radius set to " << volumeRadius_ << " mm";
00043 edm::LogInfo (MessageCategory_) << "Volume Z set to " << volumeZ_ << " mm";
00044 edm::LogInfo (MessageCategory_) << "Lower pT Threshold set to " << pTcut_1_ << " GeV";
00045 edm::LogInfo (MessageCategory_) << "Upper pT Threshold set to " << pTcut_2_ << " GeV";
00046
00047
00048 produces< std::vector<PileupSummaryInfo> >();
00049
00050 }
00051
00052
00053 void PileupInformation::produce(edm::Event &event, const edm::EventSetup & setup)
00054 {
00055
00056 std::auto_ptr<std::vector<PileupSummaryInfo> > PSIVector(new std::vector<PileupSummaryInfo>);
00057
00058 edm::Handle< PileupMixingContent > MixingPileup;
00059 event.getByLabel(PileupInfoLabel_, MixingPileup);
00060
00061 std::vector<int> BunchCrossings;
00062 std::vector<int> Interactions_Xing;
00063 std::vector<float> TrueInteractions_Xing;
00064
00065 const PileupMixingContent* MixInfo = MixingPileup.product();
00066
00067 if(MixInfo) {
00068
00069 const std::vector<int> bunchCrossing = MixInfo->getMix_bunchCrossing();
00070 const std::vector<int> interactions = MixInfo->getMix_Ninteractions();
00071 const std::vector<float> TrueInteractions = MixInfo->getMix_TrueInteractions();
00072
00073
00074 for(int ib=0; ib<(int)bunchCrossing.size(); ++ib){
00075
00076 BunchCrossings.push_back(bunchCrossing[ib]);
00077 Interactions_Xing.push_back(interactions[ib]);
00078 TrueInteractions_Xing.push_back(TrueInteractions[ib]);
00079 }
00080 }
00081 else{
00082
00083
00084 edm::Handle<CrossingFrame<SimVertex> > cfSimVertexes;
00085 if( event.getByLabel("mix", simHitLabel_, cfSimVertexes) ) {
00086
00087
00088 simVertexes_ = std::auto_ptr<MixCollection<SimVertex> >( new MixCollection<SimVertex>(cfSimVertexes.product()) );
00089
00090 int index = 0;
00091
00092
00093 EncodedEventIdToIndex vertexId;
00094 EncodedEventId oldEventId;
00095 unsigned int oldVertexId = 0;
00096 int oldBX = -1000;
00097 int oldEvent = 0;
00098
00099 std::vector<int> BunchCrossings2;
00100 std::list<int> Interactions_Xing2;
00101
00102
00103
00104 for (MixCollection<SimVertex>::MixItr iterator = simVertexes_->begin(); iterator != simVertexes_->end(); ++iterator, ++index)
00105 {
00106
00107 if (!index || iterator->eventId() != oldEventId)
00108 {
00109 if(iterator->eventId().bunchCrossing()==0 && iterator->eventId().event()==0){
00110 continue;
00111 }
00112 if(iterator->eventId().bunchCrossing() != oldBX) {
00113 BunchCrossings2.push_back(iterator->eventId().bunchCrossing());
00114 Interactions_Xing2.push_back(iterator->eventId().event());
00115 oldBX = iterator->eventId().bunchCrossing();
00116 oldEvent = iterator->eventId().event();
00117 }
00118 else { Interactions_Xing2.pop_back();
00119 Interactions_Xing2.push_back(iterator->eventId().event());
00120 oldEvent = iterator->eventId().event();
00121 }
00122
00123
00124 oldEventId = iterator->eventId();
00125 oldVertexId = iterator->vertexId();
00126 continue;
00127 }
00128
00129 }
00130
00131 std::vector<int>::iterator viter;
00132 std::list<int>::iterator liter = Interactions_Xing2.begin();
00133
00134 for(viter = BunchCrossings2.begin(); viter != BunchCrossings2.end(); ++viter, ++liter){
00135
00136 BunchCrossings.push_back((*viter));
00137 Interactions_Xing.push_back((*liter));
00138 TrueInteractions_Xing.push_back(-1.);
00139 }
00140
00141 }
00142 }
00143
00144
00145
00146
00147
00148 zpositions.clear();
00149 sumpT_lowpT.clear();
00150 sumpT_highpT.clear();
00151 ntrks_lowpT.clear();
00152 ntrks_highpT.clear();
00153 event_index_.clear();
00154
00155 int lastEvent = 0;
00156
00157 int lastBunchCrossing = 0;
00158
00159 edm::Handle<TrackingParticleCollection> mergedPH;
00160 edm::Handle<TrackingVertexCollection> mergedVH;
00161
00162 bool HaveTrackingParticles = false;
00163
00164 TrackingVertexCollection::const_iterator iVtx;
00165 TrackingVertexCollection::const_iterator iVtxTest;
00166 TrackingParticleCollection::const_iterator iTrackTest;
00167
00168
00169 if(event.getByLabel(trackingTruth_, mergedPH) && event.getByLabel(trackingTruth_, mergedVH)) {
00170
00171 HaveTrackingParticles = true;
00172
00173 iVtxTest = mergedVH->begin();
00174 iTrackTest = mergedPH->begin();
00175 }
00176
00177 int nminb_vtx = 0;
00178
00179
00180
00181 std::vector<int>::iterator BXIter;
00182 std::vector<int>::iterator InteractionsIter = Interactions_Xing.begin();
00183 std::vector<float>::iterator TInteractionsIter = TrueInteractions_Xing.begin();
00184
00185
00186
00187 for( BXIter = BunchCrossings.begin(); BXIter != BunchCrossings.end(); ++BXIter, ++InteractionsIter, ++TInteractionsIter) {
00188
00189
00190
00191 if(HaveTrackingParticles) {
00192
00193 for (iVtx = iVtxTest; iVtx != mergedVH->end(); ++iVtx) {
00194
00195 if(iVtx->eventId().bunchCrossing() == (*BXIter) ) {
00196
00197 if(iVtx->eventId().event() != lastEvent) {
00198
00199
00200
00201 float zpos = 0.;
00202 zpos = iVtx->position().z();
00203 zpositions.push_back(zpos);
00204 sumpT_lowpT.push_back(0.);
00205 sumpT_highpT.push_back(0.);
00206 ntrks_lowpT.push_back(0);
00207 ntrks_highpT.push_back(0);
00208
00209 lastEvent = iVtx->eventId().event();
00210 iVtxTest = --iVtx;
00211
00212
00213
00214 event_index_.insert(myindex::value_type(lastEvent,nminb_vtx));
00215
00216 ++nminb_vtx;
00217
00218 continue;
00219 }
00220 }
00221 }
00222
00223
00224
00225 for (TrackingParticleCollection::const_iterator iTrack = iTrackTest; iTrack != mergedPH->end(); ++iTrack)
00226 {
00227 bool FoundTrk = false;
00228
00229 float zpos=0.;
00230
00231 if(iTrack->eventId().bunchCrossing() == (*BXIter) && iTrack->eventId().event() > 0 )
00232 {
00233 FoundTrk = true;
00234 int correct_index = event_index_[iTrack->eventId().event()];
00235
00236
00237
00238 zpos = zpositions[correct_index];
00239 if(iTrack->matchedHit()>0) {
00240 if(fabs(iTrack->parentVertex()->position().z()-zpos)<0.1) {
00241
00242 float Tpx = iTrack->p4().px();
00243 float Tpy = iTrack->p4().py();
00244 float TpT = sqrt(Tpx*Tpx + Tpy*Tpy);
00245 if( TpT>pTcut_1_ ) {
00246 sumpT_lowpT[correct_index]+=TpT;
00247 ++ntrks_lowpT[correct_index];
00248 }
00249 if( TpT>pTcut_2_ ){
00250 sumpT_highpT[correct_index]+=TpT;
00251 ++ntrks_highpT[correct_index];
00252 }
00253 }
00254 }
00255 }
00256 else{
00257 if(FoundTrk) {
00258
00259 iTrackTest = --iTrack;
00260 --iTrackTest;
00261 break;
00262 }
00263
00264 }
00265
00266 }
00267
00268 }
00269
00270
00271
00272
00273
00274
00275
00276 if(!HaveTrackingParticles) {
00277
00278 zpositions.push_back(-999.);
00279 sumpT_lowpT.push_back(0.);
00280 sumpT_highpT.push_back(0.);
00281 ntrks_lowpT.push_back(0);
00282 ntrks_highpT.push_back(0);
00283
00284 }
00285
00286 PileupSummaryInfo PSI_bunch = PileupSummaryInfo(
00287 (*InteractionsIter),
00288 zpositions,
00289 sumpT_lowpT,
00290 sumpT_highpT,
00291 ntrks_lowpT,
00292 ntrks_highpT,
00293 (*BXIter),
00294 (*TInteractionsIter)
00295 );
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309 PSIVector->push_back(PSI_bunch);
00310
00311 if(HaveTrackingParticles) lastBunchCrossing = iVtx->eventId().bunchCrossing();
00312
00313 event_index_.clear();
00314 zpositions.clear();
00315 sumpT_lowpT.clear();
00316 sumpT_highpT.clear();
00317 ntrks_lowpT.clear();
00318 ntrks_highpT.clear();
00319 nminb_vtx = 0;
00320 lastEvent=0;
00321
00322
00323 }
00324
00325
00326
00327 event.put(PSIVector);
00328
00329
00330 }
00331
00332
00333 DEFINE_FWK_MODULE(PileupInformation);