Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:42

0001 // -*- C++ -*-
0002 //
0003 // Class:      TestSuite
0004 //
0005 /**\class TestSuite
0006 
0007    Description: test suite for Mixing Module
0008 
0009 */
0010 //
0011 // Original Author:  Ursula Berthon
0012 //         Created:  Fri Sep 23 11:38:38 CEST 2005
0013 //
0014 //
0015 
0016 #include "Validation/Mixing/interface/TestSuite.h"
0017 
0018 // system include files
0019 #include <memory>
0020 #include <utility>
0021 
0022 // user include files
0023 
0024 #include "DQMServices/Core/interface/DQMStore.h"
0025 #include "TFile.h"
0026 
0027 using namespace edm;
0028 
0029 TestSuite::TestSuite(const edm::ParameterSet &iConfig)
0030     : filename_(iConfig.getParameter<std::string>("fileName")),
0031       bunchcr_(iConfig.getParameter<int>("BunchNr")),
0032       minbunch_(iConfig.getParameter<int>("minBunch")),
0033       maxbunch_(iConfig.getParameter<int>("maxBunch")),
0034       dbe_(nullptr),
0035       cfTrackToken_(consumes<CrossingFrame<SimTrack>>(iConfig.getParameter<edm::InputTag>("cfTrackTag"))),
0036       cfVertexToken_(consumes<CrossingFrame<SimTrack>>(iConfig.getParameter<edm::InputTag>("cfVertexTag"))),
0037       g4SimHits_Token_(consumes<CrossingFrame<PSimHit>>(edm::InputTag("mix", "g4SimHitsTrackerHitsTECLowTof"))),
0038       g4SimHits_Ecal_Token_(consumes<CrossingFrame<PCaloHit>>(edm::InputTag("mix", "g4SimHitsEcalHitsEB"))),
0039       g4SimHits_HCal_Token_(consumes<CrossingFrame<PCaloHit>>(edm::InputTag("mix", "g4SimHitsHcalHits"))) {
0040   std::cout << "Constructed testSuite , bunchcr " << bunchcr_ << " filename: " << filename_ << std::endl;
0041 }
0042 
0043 TestSuite::~TestSuite() {
0044   // do anything here that needs to be done at desctruction time
0045   // (e.g. close files, deallocate resources etc.)
0046 }
0047 
0048 void TestSuite::beginJob() {
0049   // get hold of back-end interface
0050   dbe_ = Service<DQMStore>().operator->();
0051   dbe_->setCurrentFolder("MixingV/Mixing");
0052 }
0053 
0054 void TestSuite::endJob() {
0055   if (!filename_.empty() && dbe_)
0056     dbe_->save(filename_);
0057 }
0058 
0059 // ------------ method called to analyze the data  ------------
0060 void TestSuite::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0061   using namespace edm;
0062 
0063   // Get input
0064   edm::Handle<CrossingFrame<SimTrack>> cf_track;
0065   edm::Handle<CrossingFrame<SimVertex>> cf_vertex;
0066   edm::Handle<CrossingFrame<PSimHit>> cf_simhit;
0067   edm::Handle<CrossingFrame<PCaloHit>> cf_calohitEcal;
0068   edm::Handle<CrossingFrame<PCaloHit>> cf_calohitHcal;
0069   std::string subdetTracker("g4SimHitsTrackerHitsTECLowTof");
0070   std::string ecalsubdet("g4SimHitsEcalHitsEB");
0071   std::string hcalsubdet("g4SimHitsHcalHits");
0072   iEvent.getByToken(cfTrackToken_, cf_track);
0073   iEvent.getByToken(cfVertexToken_, cf_vertex);
0074   iEvent.getByToken(g4SimHits_Token_, cf_simhit);
0075   iEvent.getByToken(g4SimHits_Ecal_Token_, cf_calohitEcal);
0076   iEvent.getByToken(g4SimHits_HCal_Token_, cf_calohitHcal);
0077 
0078   // use MixCollection and its iterator
0079   // Please note that bunch() and getTrigger() are methods of the iterator
0080   // itself while operator-> points to the templated objects!!!!
0081 
0082   // track histo
0083   char histotracks[30], sighistotracks[30], histotracksindsig[30], histotracksind[30];
0084   sprintf(histotracks, "Tracks_bcr_%d", bunchcr_);
0085   sprintf(sighistotracks, "SignalTracks_bcr_%d", bunchcr_);
0086   sprintf(histotracksind, "VtxPointers_%d", bunchcr_);
0087   sprintf(histotracksindsig, "VtxPointers_signal_%d", bunchcr_);
0088   MonitorElement *trhist =
0089       dbe_->book1D(histotracks, "Bunchcrossings", maxbunch_ - minbunch_ + 1, minbunch_, maxbunch_ + 1);
0090   MonitorElement *trhistsig =
0091       dbe_->book1D(sighistotracks, "Bunchcrossings", maxbunch_ - minbunch_ + 1, minbunch_, maxbunch_ + 1);
0092   MonitorElement *trindhist = dbe_->book1D(histotracksind, "Track to Vertex indices", 100, 0, 500);
0093   MonitorElement *trindhistsig = dbe_->book1D(histotracksindsig, "Signal Track to Vertex indices", 100, 0, 500);
0094   std::unique_ptr<MixCollection<SimTrack>> col1(new MixCollection<SimTrack>(cf_track.product()));
0095   MixCollection<SimTrack>::iterator cfi1;
0096   for (cfi1 = col1->begin(); cfi1 != col1->end(); cfi1++) {
0097     if (cfi1.getTrigger() == 0) {
0098       trhist->Fill(cfi1.bunch());
0099       trindhist->Fill(cfi1->vertIndex());
0100     } else {
0101       trindhistsig->Fill(cfi1->vertIndex());
0102       trhistsig->Fill(cfi1.bunch());
0103     }
0104   }
0105 
0106   // vertex histo
0107   char histovertices[30], sighistovertices[30], histovertexindices[30], histovertexindicessig[30];
0108   sprintf(histovertices, "Vertices_bcr_%d", bunchcr_);
0109   sprintf(sighistovertices, "SignalVertices_bcr_%d", bunchcr_);
0110   sprintf(histovertexindices, "TrackPointers_%d", bunchcr_);
0111   sprintf(histovertexindicessig, "TrackPointers_signal_%d", bunchcr_);
0112   MonitorElement *vtxhist =
0113       dbe_->book1D(histovertices, "Bunchcrossings", maxbunch_ - minbunch_ + 1, minbunch_, maxbunch_ + 1);
0114   MonitorElement *vtxhistsig =
0115       dbe_->book1D(sighistovertices, "Bunchcrossings", maxbunch_ - minbunch_ + 1, minbunch_, maxbunch_ + 1);
0116   MonitorElement *vtxindhist = dbe_->book1D(histovertexindices, "Vertex to Track Indices", 100, 0, 300);
0117   MonitorElement *vtxindhistsig = dbe_->book1D(histovertexindicessig, "Signal Vertex to Track Indices", 100, 0, 300);
0118   std::unique_ptr<MixCollection<SimVertex>> col2(new MixCollection<SimVertex>(cf_vertex.product()));
0119   MixCollection<SimVertex>::iterator cfi2;
0120   for (cfi2 = col2->begin(); cfi2 != col2->end(); cfi2++) {
0121     if (cfi2.getTrigger() == 0) {
0122       vtxhist->Fill(cfi2.bunch());
0123       if (!cfi2->noParent())
0124         vtxindhist->Fill(cfi2->parentIndex());
0125     } else {
0126       vtxhistsig->Fill(cfi2.bunch());
0127       if (!cfi2->noParent())
0128         vtxindhistsig->Fill(cfi2->parentIndex());
0129     }
0130   }
0131 
0132   // tracker
0133   int bsp = cf_simhit->getBunchSpace();
0134   char tof[30];
0135 
0136   sprintf(tof, "TrackerHit_Tof_bcr_%d", bunchcr_);
0137   MonitorElement *tofhist =
0138       dbe_->book1D(tof, "TrackerHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0139   sprintf(tof, "SignalTrackerHit_Tof_bcr_%d", bunchcr_);
0140   MonitorElement *tofhist_sig =
0141       dbe_->book1D(tof, "TrackerHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0142   std::unique_ptr<MixCollection<PSimHit>> colsh(new MixCollection<PSimHit>(cf_simhit.product()));
0143   MixCollection<PSimHit>::iterator cfish;
0144   for (cfish = colsh->begin(); cfish != colsh->end(); cfish++) {
0145     if (cfish.getTrigger()) {
0146       tofhist_sig->Fill(cfish->timeOfFlight());
0147     } else {
0148       tofhist->Fill(cfish->timeOfFlight());
0149     }
0150   }
0151 
0152   // Ecal
0153   sprintf(tof, "EcalEBHit_Tof_bcr_%d", bunchcr_);
0154   MonitorElement *tofecalhist =
0155       dbe_->book1D(tof, "EcalEBHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0156   sprintf(tof, "SignalEcalEBHit_Tof_bcr_%d", bunchcr_);
0157   MonitorElement *tofecalhist_sig =
0158       dbe_->book1D(tof, "EcalEBHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0159   //    std::string ecalsubdet("EcalHitsEB");
0160   std::unique_ptr<MixCollection<PCaloHit>> colecal(new MixCollection<PCaloHit>(cf_calohitEcal.product()));
0161   MixCollection<PCaloHit>::iterator cfiecal;
0162   for (cfiecal = colecal->begin(); cfiecal != colecal->end(); cfiecal++) {
0163     if (cfiecal.getTrigger())
0164       tofecalhist_sig->Fill(cfiecal->time());
0165     else
0166       tofecalhist->Fill(cfiecal->time());
0167   }
0168 
0169   // Hcal
0170   sprintf(tof, "HcalHit_Tof_bcr_%d", bunchcr_);
0171   MonitorElement *tofhcalhist =
0172       dbe_->book1D(tof, "HcalHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0173   sprintf(tof, "SignalHcalHit_Tof_bcr_%d", bunchcr_);
0174   MonitorElement *tofhcalhist_sig =
0175       dbe_->book1D(tof, "HcalHit_ToF", 100, float(bsp * minbunch_), float(bsp * maxbunch_) + 50.);
0176   //    std::string hcalsubdet("HcalHits");
0177   std::unique_ptr<MixCollection<PCaloHit>> colhcal(new MixCollection<PCaloHit>(cf_calohitHcal.product()));
0178   MixCollection<PCaloHit>::iterator cfihcal;
0179 
0180   for (cfihcal = colhcal->begin(); cfihcal != colhcal->end(); cfihcal++) {
0181     if (cfihcal.getTrigger())
0182       tofhcalhist_sig->Fill(cfihcal->time());
0183     else
0184       tofhcalhist->Fill(cfihcal->time());
0185   }
0186 }