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 int oldBX = -1000;
00096
00097 std::vector<int> BunchCrossings2;
00098 std::list<int> Interactions_Xing2;
00099
00100
00101
00102 for (MixCollection<SimVertex>::MixItr iterator = simVertexes_->begin(); iterator != simVertexes_->end(); ++iterator, ++index)
00103 {
00104
00105 if (!index || iterator->eventId() != oldEventId)
00106 {
00107 if(iterator->eventId().bunchCrossing()==0 && iterator->eventId().event()==0){
00108 continue;
00109 }
00110 if(iterator->eventId().bunchCrossing() != oldBX) {
00111 BunchCrossings2.push_back(iterator->eventId().bunchCrossing());
00112 Interactions_Xing2.push_back(iterator->eventId().event());
00113 oldBX = iterator->eventId().bunchCrossing();
00114 }
00115 else { Interactions_Xing2.pop_back();
00116 Interactions_Xing2.push_back(iterator->eventId().event());
00117 }
00118
00119
00120 oldEventId = iterator->eventId();
00121 continue;
00122 }
00123
00124 }
00125
00126 std::vector<int>::iterator viter;
00127 std::list<int>::iterator liter = Interactions_Xing2.begin();
00128
00129 for(viter = BunchCrossings2.begin(); viter != BunchCrossings2.end(); ++viter, ++liter){
00130
00131 BunchCrossings.push_back((*viter));
00132 Interactions_Xing.push_back((*liter));
00133 TrueInteractions_Xing.push_back(-1.);
00134 }
00135
00136 }
00137 }
00138
00139
00140
00141
00142
00143 zpositions.clear();
00144 sumpT_lowpT.clear();
00145 sumpT_highpT.clear();
00146 ntrks_lowpT.clear();
00147 ntrks_highpT.clear();
00148 event_index_.clear();
00149
00150 int lastEvent = 0;
00151
00152
00153
00154 edm::Handle<TrackingParticleCollection> mergedPH;
00155 edm::Handle<TrackingVertexCollection> mergedVH;
00156
00157 bool HaveTrackingParticles = false;
00158
00159 TrackingVertexCollection::const_iterator iVtx;
00160 TrackingVertexCollection::const_iterator iVtxTest;
00161 TrackingParticleCollection::const_iterator iTrackTest;
00162
00163
00164 if(event.getByLabel(trackingTruth_, mergedPH) && event.getByLabel(trackingTruth_, mergedVH)) {
00165
00166 HaveTrackingParticles = true;
00167
00168 iVtxTest = mergedVH->begin();
00169 iTrackTest = mergedPH->begin();
00170 }
00171
00172 int nminb_vtx = 0;
00173
00174
00175
00176 std::vector<int>::iterator BXIter;
00177 std::vector<int>::iterator InteractionsIter = Interactions_Xing.begin();
00178 std::vector<float>::iterator TInteractionsIter = TrueInteractions_Xing.begin();
00179
00180
00181
00182 for( BXIter = BunchCrossings.begin(); BXIter != BunchCrossings.end(); ++BXIter, ++InteractionsIter, ++TInteractionsIter) {
00183
00184
00185
00186 if(HaveTrackingParticles) {
00187
00188 for (iVtx = iVtxTest; iVtx != mergedVH->end(); ++iVtx) {
00189
00190 if(iVtx->eventId().bunchCrossing() == (*BXIter) ) {
00191
00192 if(iVtx->eventId().event() != lastEvent) {
00193
00194
00195
00196 float zpos = 0.;
00197 zpos = iVtx->position().z();
00198 zpositions.push_back(zpos);
00199 sumpT_lowpT.push_back(0.);
00200 sumpT_highpT.push_back(0.);
00201 ntrks_lowpT.push_back(0);
00202 ntrks_highpT.push_back(0);
00203
00204 lastEvent = iVtx->eventId().event();
00205 iVtxTest = --iVtx;
00206
00207
00208
00209 event_index_.insert(myindex::value_type(lastEvent,nminb_vtx));
00210
00211 ++nminb_vtx;
00212
00213 continue;
00214 }
00215 }
00216 }
00217
00218
00219
00220 for (TrackingParticleCollection::const_iterator iTrack = iTrackTest; iTrack != mergedPH->end(); ++iTrack)
00221 {
00222 bool FoundTrk = false;
00223
00224 float zpos=0.;
00225
00226 if(iTrack->eventId().bunchCrossing() == (*BXIter) && iTrack->eventId().event() > 0 )
00227 {
00228 FoundTrk = true;
00229 int correct_index = event_index_[iTrack->eventId().event()];
00230
00231
00232
00233 zpos = zpositions[correct_index];
00234 if(iTrack->matchedHit()>0) {
00235 if(fabs(iTrack->parentVertex()->position().z()-zpos)<0.1) {
00236
00237 float Tpx = iTrack->p4().px();
00238 float Tpy = iTrack->p4().py();
00239 float TpT = sqrt(Tpx*Tpx + Tpy*Tpy);
00240 if( TpT>pTcut_1_ ) {
00241 sumpT_lowpT[correct_index]+=TpT;
00242 ++ntrks_lowpT[correct_index];
00243 }
00244 if( TpT>pTcut_2_ ){
00245 sumpT_highpT[correct_index]+=TpT;
00246 ++ntrks_highpT[correct_index];
00247 }
00248 }
00249 }
00250 }
00251 else{
00252 if(FoundTrk) {
00253
00254 iTrackTest = --iTrack;
00255 --iTrackTest;
00256 break;
00257 }
00258
00259 }
00260
00261 }
00262
00263 }
00264
00265
00266
00267
00268
00269
00270
00271 if(!HaveTrackingParticles) {
00272
00273 zpositions.push_back(-999.);
00274 sumpT_lowpT.push_back(0.);
00275 sumpT_highpT.push_back(0.);
00276 ntrks_lowpT.push_back(0);
00277 ntrks_highpT.push_back(0);
00278
00279 }
00280
00281 PileupSummaryInfo PSI_bunch = PileupSummaryInfo(
00282 (*InteractionsIter),
00283 zpositions,
00284 sumpT_lowpT,
00285 sumpT_highpT,
00286 ntrks_lowpT,
00287 ntrks_highpT,
00288 (*BXIter),
00289 (*TInteractionsIter)
00290 );
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304 PSIVector->push_back(PSI_bunch);
00305
00306
00307
00308 event_index_.clear();
00309 zpositions.clear();
00310 sumpT_lowpT.clear();
00311 sumpT_highpT.clear();
00312 ntrks_lowpT.clear();
00313 ntrks_highpT.clear();
00314 nminb_vtx = 0;
00315 lastEvent=0;
00316
00317
00318 }
00319
00320
00321
00322 event.put(PSIVector);
00323
00324
00325 }
00326
00327
00328 DEFINE_FWK_MODULE(PileupInformation);