Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:37

0001 // -*- C++ -*-
0002 //
0003 // Class:      TestMix
0004 //
0005 /**\class TestMix
0006 
0007  Description: test of Mixing Module
0008 
0009 */
0010 //
0011 // Original Author:  Ursula Berthon
0012 //         Created:  Fri Sep 23 11:38:38 CEST 2005
0013 //
0014 //
0015 
0016 // system include files
0017 #include <memory>
0018 #include <utility>
0019 #include <iostream>
0020 
0021 // user include files
0022 #include "FWCore/Framework/interface/Frameworkfwd.h"
0023 
0024 #include "FWCore/Framework/interface/Event.h"
0025 #include "FWCore/Framework/interface/MakerMacros.h"
0026 
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 
0029 #include "TestMix.h"
0030 
0031 using namespace edm;
0032 
0033 TestMix::TestMix(const edm::ParameterSet& iConfig) : level_(iConfig.getUntrackedParameter<int>("PrintLevel")) {
0034   std::cout << "Constructed testMix , level " << level_ << std::endl;
0035 
0036   track_containers_.push_back("g4SimHitsTrackerHitsTECHighTof");
0037   track_containers_.push_back("g4SimHitsTrackerHitsTECLowTof");
0038 
0039   track_containers2_.push_back("g4SimHitsTrackerHitsTECLowTof");
0040   track_containers2_.push_back("g4SimHitsTrackerHitsTECHighTof");
0041 
0042   edm::InputTag tag = edm::InputTag("mix", "g4SimHits");
0043 
0044   SimTrackToken_ = consumes<CrossingFrame<SimTrack>>(tag);
0045   SimVertexToken_ = consumes<CrossingFrame<SimVertex>>(tag);
0046 
0047   tag = edm::InputTag("mix", "g4SimHitsTrackerHitsTECHighTof");
0048   TrackerToken0_ = consumes<CrossingFrame<PSimHit>>(tag);
0049 
0050   tag = edm::InputTag("mix", "g4SimHitsEcalHitsEB");
0051   CaloToken1_ = consumes<CrossingFrame<PCaloHit>>(tag);
0052 
0053   tag = edm::InputTag("mix", track_containers_[0]);
0054   TrackerToken1_ = consumes<CrossingFrame<PSimHit>>(tag);
0055 
0056   tag = edm::InputTag("mix", track_containers_[1]);
0057   TrackerToken2_ = consumes<CrossingFrame<PSimHit>>(tag);
0058 
0059   tag = edm::InputTag("mix", track_containers2_[0]);
0060   TrackerToken3_ = consumes<CrossingFrame<PSimHit>>(tag);
0061 
0062   tag = edm::InputTag("mix", track_containers2_[1]);
0063   TrackerToken4_ = consumes<CrossingFrame<PSimHit>>(tag);
0064 
0065   tag = edm::InputTag("mix", "generatorSmeared");
0066   HepMCToken_ = consumes<CrossingFrame<HepMCProduct>>(tag);
0067 }
0068 
0069 TestMix::~TestMix() {
0070   // do anything here that needs to be done at desctruction time
0071   // (e.g. close files, deallocate resources etc.)
0072 }
0073 
0074 //
0075 // member functions
0076 //
0077 
0078 // ------------ method called to analyze the data  ------------
0079 void TestMix::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0080   using namespace edm;
0081   bool got;
0082   int count = 0;
0083 
0084   // test accesses to CrossingFrame
0085   // attention: operator-> returns the templated object, but
0086   // bunch() and getTrigger() are methods of the iterator itself!
0087 
0088   // test access to SimTracks directly in CrossingFrame
0089 
0090   edm::Handle<CrossingFrame<SimTrack>> cf_simtrack;
0091   bool gotTracks = iEvent.getByToken(SimTrackToken_, cf_simtrack);
0092   if (!gotTracks)
0093     std::cout << " Could not read SimTracks!!!!" << std::endl;
0094 
0095   // not pointer compatible!!!!
0096   //   if (gotTracks) {
0097   //     std::vector<SimTrack>::const_iterator first;
0098   //     std::vector<SimTrack>::const_iterator last;
0099   //     cf_simtrack->getPileups(first,last);
0100   //     unsigned int ic=0;
0101   //     for (std::vector<SimTrack>::const_iterator it=first;it!=last;it++) {
0102   //       std::cout<<" pileup SimTrack "<<ic<<" has genpart index  "<<(*it).genpartIndex()<<" vertex Index "<<(*it).vertIndex()  <<" bunchcrossing "<<cf_simtrack->getBunchCrossing(ic)<<std::endl;
0103   //       ic++;
0104   //     }
0105   //   }
0106 
0107   const std::string subdet("g4SimHitsTrackerHitsTECHighTof");
0108   edm::Handle<CrossingFrame<PSimHit>> cf_simhit;
0109   got = iEvent.getByToken(TrackerToken0_, cf_simhit);
0110   if (!got)
0111     std::cout << " Could not read SimHits with label " << subdet << "!!!!" << std::endl;
0112   else {
0113     std::cout << "\n\n=================== Starting SimHit access, subdet " << subdet
0114               << "  ===================" << std::endl;
0115 
0116     std::unique_ptr<MixCollection<PSimHit>> col(
0117         new MixCollection<PSimHit>(cf_simhit.product(), std::pair<int, int>(-1, 1)));
0118     std::cout << *(col.get()) << std::endl;
0119     MixCollection<PSimHit>::iterator cfi;
0120     for (cfi = col->begin(); cfi != col->end(); cfi++) {
0121       std::cout << " Hit " << count << " has tof " << cfi->timeOfFlight() << " trackid " << cfi->trackId()
0122                 << " bunchcr " << cfi.bunch() << " trigger " << cfi.getTrigger()
0123                 << ", from EncodedEventId: " << cfi->eventId().bunchCrossing() << " " << cfi->eventId().event()
0124                 << " bcr from MixCol " << cfi.bunch() << std::endl;
0125       //      std::cout<<" Hit: "<<(*cfi)<<std::endl;
0126       count++;
0127     }
0128   }
0129 
0130   // test access to CaloHits
0131   const std::string subdetcalo("g4SimHitsEcalHitsEB");
0132   edm::Handle<CrossingFrame<PCaloHit>> cf_calo;
0133   got = iEvent.getByToken(CaloToken1_, cf_calo);
0134   if (!got)
0135     std::cout << " Could not read CaloHits with label " << subdetcalo << "!!!!" << std::endl;
0136   else {
0137     std::cout << "\n\n=================== Starting CaloHit access, subdet " << subdetcalo
0138               << "  ===================" << std::endl;
0139     std::unique_ptr<MixCollection<PCaloHit>> colcalo(
0140         new MixCollection<PCaloHit>(cf_calo.product(), std::pair<int, int>(-1, 1)));
0141     std::cout << *(colcalo.get()) << std::endl;
0142     MixCollection<PCaloHit>::iterator cficalo;
0143     count = 0;
0144     for (cficalo = colcalo->begin(); cficalo != colcalo->end(); cficalo++) {
0145       std::cout << " CaloHit " << count << " has tof " << cficalo->time() << " trackid " << cficalo->geantTrackId()
0146                 << " bunchcr " << cficalo.bunch() << " trigger " << cficalo.getTrigger()
0147                 << ", from EncodedEventId: " << cficalo->eventId().bunchCrossing() << " " << cficalo->eventId().event()
0148                 << std::endl;
0149       //      std::cout<<" Calo Hit: "<<(*cficalo)<<std::endl;
0150       count++;
0151     }
0152   }
0153 
0154   // test access to SimTracks
0155   if (gotTracks) {
0156     std::cout << "\n=================== Starting SimTrack access ===================" << std::endl;
0157     //   edm::Handle<CrossingFrame<SimTrack> > cf_simtrack;
0158     //   iEvent.getByLabel("mix",cf_simtrack);
0159     std::unique_ptr<MixCollection<SimTrack>> col2(new MixCollection<SimTrack>(cf_simtrack.product()));
0160     MixCollection<SimTrack>::iterator cfi2;
0161     int count2 = 0;
0162     std::cout << " \nWe got " << col2->sizeSignal() << " signal tracks and " << col2->sizePileup()
0163               << " pileup tracks, total: " << col2->size() << std::endl;
0164     for (cfi2 = col2->begin(); cfi2 != col2->end(); cfi2++) {
0165       std::cout << " SimTrack " << count2 << " has genpart index  " << cfi2->genpartIndex() << " vertex Index "
0166                 << cfi2->vertIndex() << " bunchcr " << cfi2.bunch() << " trigger " << cfi2.getTrigger()
0167                 << ", from EncodedEventId: " << cfi2->eventId().bunchCrossing() << " " << cfi2->eventId().event()
0168                 << std::endl;
0169       count2++;
0170     }
0171   }
0172 
0173   // test access to SimVertices
0174   edm::Handle<CrossingFrame<SimVertex>> cf_simvtx;
0175   got = iEvent.getByToken(SimVertexToken_, cf_simvtx);
0176   if (!got)
0177     std::cout << " Could not read Simvertices !!!!" << std::endl;
0178   else {
0179     std::cout << "\n=================== Starting SimVertex access ===================" << std::endl;
0180     std::unique_ptr<MixCollection<SimVertex>> col3(new MixCollection<SimVertex>(cf_simvtx.product()));
0181     MixCollection<SimVertex>::iterator cfi3;
0182     int count3 = 0;
0183     std::cout << " \nWe got " << col3->sizeSignal() << " signal vertices and " << col3->sizePileup()
0184               << " pileup vertices, total: " << col3->size() << std::endl;
0185     for (cfi3 = col3->begin(); cfi3 != col3->end(); cfi3++) {
0186       std::cout << " SimVertex " << count3 << " has parent index  " << cfi3->parentIndex() << " bunchcr "
0187                 << cfi3.bunch() << " trigger " << cfi3.getTrigger()
0188                 << ", from EncodedEventId: " << cfi3->eventId().bunchCrossing() << " " << cfi3->eventId().event()
0189                 << std::endl;
0190       const SimVertex& myvtx = (*cfi3);
0191       std::cout << "Same with op*: " << count3 << " has parent index  " << myvtx.parentIndex() << " bunchcr "
0192                 << cfi3.bunch() << " trigger " << cfi3.getTrigger()
0193                 << ", from EncodedEventId: " << myvtx.eventId().bunchCrossing() << " " << myvtx.eventId().event()
0194                 << std::endl;
0195       count3++;
0196     }
0197   }
0198 
0199   //test MixCollection constructor with several subdetector names
0200   bool got1, got2 = false;
0201   std::unique_ptr<MixCollection<PSimHit>> all_trackhits;
0202   std::unique_ptr<MixCollection<PSimHit>> all_trackhits2;
0203   std::cout << "\n=================== Starting test for coll of several ROU-s ===================" << std::endl;
0204   //  edm::Handle<CrossingFrame<PSimHit> > cf_simhit;
0205   std::vector<const CrossingFrame<PSimHit>*> cfvec;
0206   got1 = iEvent.getByToken(TrackerToken1_, cf_simhit);
0207   if (!got1)
0208     std::cout << " Could not read SimHits with label " << track_containers_[0] << "!!!!" << std::endl;
0209   else {
0210     std::cout << "\n=================== Starting test for coll of several ROU-s ===================" << std::endl;
0211     cfvec.push_back(cf_simhit.product());
0212     std::cout << " \nFirst container " << track_containers_[0] << " Nr signals " << cf_simhit->getNrSignals()
0213               << ", Nr pileups " << cf_simhit->getNrPileups() << std::endl;
0214     got2 = iEvent.getByToken(TrackerToken2_, cf_simhit);
0215     if (got2) {
0216       cfvec.push_back(cf_simhit.product());
0217       std::cout << " \nSecond container " << track_containers_[1] << " Nr signals " << cf_simhit->getNrSignals()
0218                 << ", Nr pileups " << cf_simhit->getNrPileups() << std::endl;
0219       all_trackhits = std::make_unique<MixCollection<PSimHit>>(cfvec);
0220 
0221       std::cout << " \nFor all containers we got " << all_trackhits->sizeSignal() << " signal hits and "
0222                 << all_trackhits->sizePileup() << " pileup hits, total: " << all_trackhits->size() << std::endl;
0223 
0224       MixCollection<PSimHit>::iterator it;
0225       int ii = 0;
0226       for (it = all_trackhits->begin(); it != all_trackhits->end(); it++) {
0227         std::cout << " Hit " << ii << " of all hits has tof " << it->timeOfFlight() << " trackid " << it->trackId()
0228                   << " bunchcr " << it.bunch() << " trigger " << it.getTrigger()
0229                   << ", from EncodedEventId: " << it->eventId().bunchCrossing() << " " << it->eventId().event()
0230                   << std::endl;
0231         ii++;
0232       }
0233     }
0234   }
0235 
0236   //test the same in different order: should be the same sizes, different order
0237   MixCollection<PSimHit>::iterator it2;
0238   int ii2 = 0;
0239   std::vector<const CrossingFrame<PSimHit>*> cfvec2;
0240   got = iEvent.getByToken(TrackerToken3_, cf_simhit);
0241   if (!got)
0242     std::cout << " Could not read SimHits with label " << track_containers2_[0] << "!!!!" << std::endl;
0243   else {
0244     cfvec2.push_back(cf_simhit.product());
0245     got2 = iEvent.getByToken(TrackerToken4_, cf_simhit);
0246     if (got2) {
0247       cfvec2.push_back(cf_simhit.product());
0248       all_trackhits2 = std::make_unique<MixCollection<PSimHit>>(cfvec2);
0249       std::cout << " \nSame containers, different order: we got " << all_trackhits2->sizeSignal() << " signal hits and "
0250                 << all_trackhits2->sizePileup() << " pileup hits, total: " << all_trackhits2->size() << std::endl;
0251       for (it2 = all_trackhits2->begin(); it2 != all_trackhits2->end(); it2++) {
0252         std::cout << " Hit " << ii2 << " of all hits has tof " << it2->timeOfFlight() << " trackid " << it2->trackId()
0253                   << " bunchcr " << it2.bunch() << " trigger " << it2.getTrigger()
0254                   << ", bcr from Id: " << it2->eventId().bunchCrossing() << " evtnr in id " << it2->eventId().event()
0255                   << std::endl;
0256         ii2++;
0257       }
0258     }
0259   }
0260 
0261   //test MixCollection for HepMCProducts
0262   //at the same time test getObject method
0263   //we should have each line twice
0264   //------------------------------------
0265   edm::Handle<CrossingFrame<edm::HepMCProduct>> cf_hepmc;
0266   got = iEvent.getByToken(HepMCToken_, cf_hepmc);
0267   if (!got)
0268     std::cout << " Could not read HepMCProducts!!!!" << std::endl;
0269   else {
0270     std::unique_ptr<MixCollection<edm::HepMCProduct>> colhepmc(
0271         new MixCollection<edm::HepMCProduct>(cf_hepmc.product()));
0272     MixCollection<edm::HepMCProduct>::iterator cfihepmc;
0273     int counthepmc = 0;
0274     std::cout << " \nWe got " << colhepmc->sizeSignal() << " signal hepmc products and " << colhepmc->sizePileup()
0275               << " pileup hepmcs, total: " << colhepmc->size() << std::endl;
0276     for (cfihepmc = colhepmc->begin(); cfihepmc != colhepmc->end(); cfihepmc++) {
0277       std::cout << " edm::HepMCProduct " << counthepmc << " has event number " << cfihepmc->GetEvent()->event_number()
0278                 << ", " << cfihepmc->GetEvent()->particles_size() << " particles and "
0279                 << cfihepmc->GetEvent()->vertices_size() << " vertices,  bunchcr= " << cfihepmc.bunch()
0280                 << " trigger= " << cfihepmc.getTrigger() << " sourcetype= " << cfihepmc.getSourceType() << std::endl;
0281       HepMCProduct myprod = colhepmc->getObject(counthepmc);
0282       std::cout << "same with getObject:hepmc product   " << counthepmc << " has event number "
0283                 << myprod.GetEvent()->event_number() << ", " << myprod.GetEvent()->particles_size() << " particles and "
0284                 << myprod.GetEvent()->vertices_size() << " vertices" << std::endl;
0285       counthepmc++;
0286     }
0287   }
0288   //----------------------------------------------------------------------------
0289   //  testing special situations
0290   //----------------------------------------------------------------------------
0291 
0292   if (got2) {
0293     // test reusage of the same iterator
0294     int ii3 = 0;
0295     for (it2 = all_trackhits2->begin(); it2 != all_trackhits2->end(); it2++)
0296       ii3++;
0297     if (ii3 != ii2)
0298       std::cout << " Problem when re-using iterator!!" << std::endl;
0299     else
0300       std::cout << " \nNo problem when re-using iterator." << std::endl;
0301   }
0302   // test access to non-filled collections
0303   //cases:   0) ok, collection has members
0304   //         1) bunchrange given outside of existent bunchcrossing numbers ==>exc
0305 
0306   std::cout << "\n=================== Starting tests for abnormal conditions ===================" << std::endl;
0307 
0308   // test case 0
0309   if (got1) {
0310     std::cout << "\n[ Testing abnormal conditions case 0]Should be all ok: registry: " << all_trackhits->inRegistry()
0311               << " size: " << all_trackhits->size() << std::endl;
0312 
0313     // test case 1
0314     std::cout << "\n[ Testing abnormal conditions case 1] Should throw an exception " << std::endl;
0315     MixCollection<PSimHit>* col21 = nullptr;
0316     col21 = new MixCollection<PSimHit>(cf_simhit.product(), std::pair<int, int>(-10, 20));
0317     delete col21;
0318   }
0319 }