00001
00002 #include <memory>
00003 #include <utility>
00004
00005
00006 #include "FWCore/Framework/interface/Frameworkfwd.h"
00007 #include "FWCore/Framework/interface/EDAnalyzer.h"
00008
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/MakerMacros.h"
00011
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013
00014 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00015 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00016
00017 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00018 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00019
00020 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00021 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00022
00023 #include "Validation/Mixing/interface/MixCollectionValidation.h"
00024 #include "TFile.h"
00025 #include "DQMServices/Core/interface/DQMStore.h"
00026 #include "DQMServices/Core/interface/MonitorElement.h"
00027
00028 using namespace edm;
00029
00030 MixCollectionValidation::MixCollectionValidation(const edm::ParameterSet& iConfig): outputFile_(iConfig.getParameter<std::string>("outputFile")),
00031 minbunch_(iConfig.getParameter<int>("minBunch")),
00032 maxbunch_(iConfig.getParameter<int>("maxBunch")),
00033 verbose_(iConfig.getUntrackedParameter<bool>("verbose",false)),
00034 dbe_(0),nbin_(maxbunch_-minbunch_+1)
00035 {
00036
00037 if ( outputFile_.size() != 0 ) {
00038 edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will be saved to " << outputFile_.c_str();
00039 } else {
00040 edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will NOT be saved";
00041 }
00042
00043
00044 dbe_ = edm::Service<DQMStore>().operator->();
00045 if ( dbe_ ) {
00046 if ( verbose_ ) { dbe_->setVerbose(1); }
00047 else { dbe_->setVerbose(0); }
00048 }
00049
00050 if ( dbe_ ) {
00051 if ( verbose_ ) dbe_->showDirStructure();
00052 }
00053
00054
00055 dbe_ = Service<DQMStore>().operator->();
00056
00057 dbe_->setCurrentFolder("MixingV/Mixing");
00058
00059
00060
00061 ParameterSet ps=iConfig.getParameter<ParameterSet>("mixObjects");
00062 names_ = ps.getParameterNames();
00063
00064 for (std::vector<std::string>::iterator it=names_.begin();it!= names_.end();++it)
00065 {
00066 ParameterSet pset=ps.getParameter<ParameterSet>((*it));
00067 if (!pset.exists("type")) continue;
00068 std::string object = pset.getParameter<std::string>("type");
00069 std::vector<InputTag> tags=pset.getParameter<std::vector<InputTag> >("input");
00070
00071 if ( object == "HepMCProduct" ) {
00072
00073 std::string title = "Log10 Number of GenParticle in " + object;
00074 std::string name = "NumberOf" + object;
00075 nrHepMCProductH_ = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,0.,40.);
00076
00077 HepMCProductTags_ = tags;
00078
00079 }
00080 else if ( object == "SimTrack" ) {
00081
00082 std::string title = "Log10 Number of " + object;
00083 std::string name = "NumberOf" + object;
00084 nrSimTrackH_ = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,0.,40.);
00085
00086 SimTrackTags_ = tags;
00087
00088 }
00089 else if ( object == "SimVertex" ) {
00090
00091 std::string title = "Log10 Number of " + object;
00092 std::string name = "NumberOf" + object;
00093 nrSimVertexH_ = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,0.,40.);
00094
00095 SimVertexTags_ = tags;
00096
00097 }
00098 else if ( object == "PSimHit" ) {
00099 std::vector<std::string> subdets=pset.getParameter<std::vector<std::string> >("subdets");
00100 for (unsigned int ii=0;ii<subdets.size();ii++) {
00101
00102 std::string title = "Log10 Number of " + subdets[ii];
00103 std::string name = "NumberOf" + subdets[ii];
00104 SimHitNrmap_[subdets[ii]] = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,0.,40.);
00105
00106 title = "Time of " + subdets[ii];
00107 name = "TimeOf" + subdets[ii];
00108 SimHitTimemap_[subdets[ii]] = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,-125.,375.);
00109
00110 }
00111
00112 PSimHitTags_ = tags;
00113
00114 }
00115 else if ( object == "PCaloHit" ) {
00116 std::vector<std::string> subdets=pset.getParameter<std::vector<std::string> >("subdets");
00117 for (unsigned int ii=0;ii<subdets.size();ii++) {
00118
00119 std::string title = "Log10 Number of " + subdets[ii];
00120 std::string name = "NumberOf" + subdets[ii];
00121 CaloHitNrmap_[subdets[ii]] = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,0.,40.);
00122
00123 title = "Time of " + subdets[ii];
00124 name = "TimeOf" + subdets[ii];
00125 CaloHitTimemap_[subdets[ii]] = dbe_->bookProfile(name,title,nbin_,minbunch_,maxbunch_+1,40,-125.,375.);
00126
00127 }
00128
00129 PCaloHitTags_ = tags;
00130 }
00131 }
00132
00133 }
00134
00135 MixCollectionValidation::~MixCollectionValidation()
00136 {
00137
00138
00139
00140 }
00141
00142 void MixCollectionValidation::beginJob() {
00143 }
00144
00145
00146 void MixCollectionValidation::endJob() {
00147 if (outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00148 }
00149
00150
00151
00152
00153
00154
00155 void
00156 MixCollectionValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iConfig)
00157 {
00158
00159 using namespace edm;
00160
00161 if ( HepMCProductTags_.size() > 0 ) {
00162 bool gotHepMCProduct;
00163 edm::Handle<CrossingFrame<HepMCProduct> > crossingFrame;
00164 std::string HepMCProductLabel = HepMCProductTags_[0].label();
00165 gotHepMCProduct = iEvent.getByLabel("mix",HepMCProductLabel,crossingFrame);
00166
00167 if (gotHepMCProduct){
00168 std::auto_ptr<MixCollection<HepMCProduct> >
00169 hepMCProduct (new MixCollection<HepMCProduct>(crossingFrame.product ()));
00170 MixCollection<HepMCProduct>::MixItr hitItr;
00171
00172 fillGenParticleMulti(hitItr, hepMCProduct, nrHepMCProductH_);
00173 }
00174 }
00175
00176 if ( SimTrackTags_.size() > 0 ) {
00177 bool gotSimTrack;
00178 edm::Handle<CrossingFrame<SimTrack> > crossingFrame;
00179 std::string SimTrackLabel = SimTrackTags_[0].label();
00180 gotSimTrack = iEvent.getByLabel("mix",SimTrackLabel,crossingFrame);
00181
00182 if (gotSimTrack){
00183 std::auto_ptr<MixCollection<SimTrack> >
00184 simTracks (new MixCollection<SimTrack>(crossingFrame.product ()));
00185 MixCollection<SimTrack>::MixItr hitItr;
00186
00187 fillMultiplicity(hitItr, simTracks, nrSimTrackH_);
00188 }
00189 }
00190
00191 if ( SimVertexTags_.size() > 0 ) {
00192 bool gotSimVertex;
00193 edm::Handle<CrossingFrame<SimVertex> > crossingFrame;
00194 std::string SimVertexLabel = SimVertexTags_[0].label();
00195 gotSimVertex = iEvent.getByLabel("mix",SimVertexLabel,crossingFrame);
00196
00197 if (gotSimVertex){
00198 std::auto_ptr<MixCollection<SimVertex> >
00199 simVerteces (new MixCollection<SimVertex>(crossingFrame.product ()));
00200 MixCollection<SimVertex>::MixItr hitItr;
00201
00202 fillMultiplicity(hitItr, simVerteces, nrSimVertexH_);
00203 }
00204 }
00205
00206 if ( PSimHitTags_.size() > 0 ) {
00207
00208 edm::Handle<CrossingFrame<PSimHit> > crossingFrame;
00209
00210 for ( int i = 0; i < (int)PSimHitTags_.size(); i++ ) {
00211 bool gotPSimHit;
00212 std::string PSimHitLabel = PSimHitTags_[i].label()+PSimHitTags_[i].instance();
00213 gotPSimHit = iEvent.getByLabel("mix",PSimHitLabel,crossingFrame);
00214
00215 if (gotPSimHit){
00216 std::auto_ptr<MixCollection<PSimHit> >
00217 simHits (new MixCollection<PSimHit>(crossingFrame.product ()));
00218
00219 MixCollection<PSimHit>::MixItr hitItr;
00220
00221 fillMultiplicity(hitItr, simHits, SimHitNrmap_[PSimHitTags_[i].instance()]);
00222
00223 fillSimHitTime(hitItr, simHits, SimHitTimemap_[PSimHitTags_[i].instance()]);
00224 }
00225 }
00226 }
00227
00228 if ( PCaloHitTags_.size() > 0 ) {
00229
00230 edm::Handle<CrossingFrame<PCaloHit> > crossingFrame;
00231
00232 for ( int i = 0; i < (int)PCaloHitTags_.size(); i++ ) {
00233 bool gotPCaloHit;
00234 std::string PCaloHitLabel = PCaloHitTags_[i].label()+PCaloHitTags_[i].instance();
00235 gotPCaloHit = iEvent.getByLabel("mix",PCaloHitLabel,crossingFrame);
00236
00237 if (gotPCaloHit){
00238 std::auto_ptr<MixCollection<PCaloHit> >
00239 caloHits (new MixCollection<PCaloHit>(crossingFrame.product ()));
00240
00241 MixCollection<PCaloHit>::MixItr hitItr;
00242
00243 fillMultiplicity(hitItr, caloHits, CaloHitNrmap_[PCaloHitTags_[i].instance()]);
00244
00245 fillCaloHitTime(hitItr, caloHits, CaloHitTimemap_[PCaloHitTags_[i].instance()]);
00246 }
00247 }
00248 }
00249
00250 }
00251
00252 template<class T1, class T2> void MixCollectionValidation::fillMultiplicity(T1 & theItr_, T2 & theColl_, MonitorElement * theProfile_) {
00253
00254 std::vector<int> theMult(nbin_);
00255
00256 for ( theItr_ = theColl_->begin() ; theItr_ != theColl_->end() ; ++theItr_) {
00257
00258 int bunch = (*theItr_).eventId().bunchCrossing();
00259 int index = bunch - minbunch_;
00260 if ( index >= 0 && index < nbin_ ) { theMult[index] += 1; }
00261 else { edm::LogWarning("MixCollectionValidation") << "fillMultiplicity: bunch number " << bunch << " out of range"; }
00262
00263 }
00264
00265 for ( int i = 0; i < nbin_; i++ ) {
00266 theProfile_->Fill(float(i+minbunch_+0.5),std::log10(std::max(float(0.1),float(theMult[i]))));
00267 }
00268 }
00269
00270
00271 template<class T1, class T2> void MixCollectionValidation::fillGenParticleMulti(T1 & theItr_, T2 & theColl_, MonitorElement * theProfile_) {
00272
00273 std::vector<int> theMult(nbin_);
00274
00275 for ( theItr_ = theColl_->begin() ; theItr_ != theColl_->end() ; ++theItr_) {
00276
00277 int bunch = theItr_.bunch();
00278 int index = bunch - minbunch_;
00279 if ( index >= 0 && index < nbin_ ) { theMult[index] += (*theItr_).GetEvent()->particles_size(); }
00280 else { edm::LogWarning("MixCollectionValidation") << "fillMultiplicity: bunch number " << bunch << " out of range"; }
00281
00282 }
00283
00284 for ( int i = 0; i < nbin_; i++ ) {
00285 theProfile_->Fill(float(i+minbunch_+0.5),std::log10(std::max(float(0.1),float(theMult[i]))));
00286 }
00287 }
00288
00289 template<class T1, class T2> void MixCollectionValidation::fillSimHitTime(T1 & theItr_, T2 & theColl_, MonitorElement * theProfile_) {
00290
00291 for ( theItr_ = theColl_->begin() ; theItr_ != theColl_->end() ; ++theItr_) {
00292
00293 int bunch = (*theItr_).eventId().bunchCrossing();
00294 float time = (*theItr_).timeOfFlight();
00295 int index = bunch - minbunch_;
00296 if ( index >= 0 && index < nbin_ ) { theProfile_->Fill(float(bunch+0.5),time); }
00297 else { edm::LogWarning("MixCollectionValidation") << "fillSimHitTime: bunch number " << bunch << " out of range"; }
00298
00299 }
00300
00301 }
00302
00303 template<class T1, class T2> void MixCollectionValidation::fillCaloHitTime(T1 & theItr_, T2 & theColl_, MonitorElement * theProfile_) {
00304
00305 for ( theItr_ = theColl_->begin() ; theItr_ != theColl_->end() ; ++theItr_) {
00306
00307 int bunch = (*theItr_).eventId().bunchCrossing();
00308 float time = (*theItr_).time();
00309 int index = bunch - minbunch_;
00310 if ( index >= 0 && index < nbin_ ) { theProfile_->Fill(float(bunch+0.5),time); }
00311 else { edm::LogWarning("MixCollectionValidation") << "fillCaloHitTime: bunch number " << bunch << " out of range"; }
00312
00313 }
00314
00315 }