Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*class BasicHepMCHeavyIonValidation
0002  *  
0003  *  Class to fill dqm monitor elements from existing EDM file
0004  *  Quan Wang - 04/2013
0005  */
0006 
0007 #include "Validation/EventGenerator/interface/BasicHepMCHeavyIonValidation.h"
0008 
0009 #include "CLHEP/Units/defs.h"
0010 #include "CLHEP/Units/PhysicalConstants.h"
0011 #include "Validation/EventGenerator/interface/DQMHelper.h"
0012 
0013 using namespace edm;
0014 
0015 BasicHepMCHeavyIonValidation::BasicHepMCHeavyIonValidation(const edm::ParameterSet& iPSet)
0016     : wmanager_(iPSet, consumesCollector()), hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")) {
0017   QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug", false);
0018 
0019   hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
0020 }
0021 
0022 BasicHepMCHeavyIonValidation::~BasicHepMCHeavyIonValidation() {}
0023 
0024 void BasicHepMCHeavyIonValidation::bookHistograms(DQMStore::IBooker& i, edm::Run const&, edm::EventSetup const&) {
0025   ///Setting the DQM top directories
0026   DQMHelper dqm(&i);
0027   i.setCurrentFolder("Generator/HeavyIon");
0028 
0029   // Number of analyzed events
0030   nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "", "Number of Events");
0031 
0032   ///Booking the ME's
0033   Ncoll_hard =
0034       dqm.book1dHisto("Ncoll_hard", "Ncoll_hard", 700, 0, 700, "Number of hard scatterings", "Number of Events");
0035   Npart_proj =
0036       dqm.book1dHisto("Npart_proj", "Npart_proj", 250, 0, 250, "Number of projectile participants", "Number of Events");
0037   Npart_targ =
0038       dqm.book1dHisto("Npart_targ", "Npart_targ", 250, 0, 250, "Number of target participants", "Number of Events");
0039   Ncoll = dqm.book1dHisto("Ncoll", "Ncoll", 700, 0, 700, "Number of N-N collisions", "Number of Events");
0040   N_Nwounded_collisions = dqm.book1dHisto("N_Nwounded_collisions",
0041                                           "N_Nwounded_collisions",
0042                                           250,
0043                                           0,
0044                                           250,
0045                                           "Number of N-N wounded collisions",
0046                                           "Number of Events");
0047   Nwounded_N_collisions = dqm.book1dHisto("Nwounded_N_collisions",
0048                                           "Nwounded_N_collisions",
0049                                           250,
0050                                           0,
0051                                           250,
0052                                           "Number of N wounded-N collisions",
0053                                           "Number of Events");
0054   Nwounded_Nwounded_collisions = dqm.book1dHisto("Nwounded_Nwounded_collisions",
0055                                                  "Nwounded_Nwounded_collisions",
0056                                                  250,
0057                                                  0,
0058                                                  250,
0059                                                  "Number of N wounded-N wounded collisions",
0060                                                  "Number of Events");
0061   spectator_neutrons = dqm.book1dHisto(
0062       "spectator_neutrons", "spectator_neutrons", 250, 0, 250, "Number of spectator neutrons", "Number of Events");
0063   spectator_protons = dqm.book1dHisto(
0064       "spectator_protons", "spectator_protons", 250, 0, 250, "Number of spectator protons", "Number of Events");
0065   impact_parameter = dqm.book1dHisto(
0066       "impact_parameter", "impact_parameter", 50, 0, 50, "Empact parameter of collision (fm)", "Number of Events");
0067   event_plane_angle = dqm.book1dHisto(
0068       "event_plane_angle", "event_plane_angle", 200, 0, 2 * CLHEP::pi, "#phi_{event plane} (rad)", "Number of Events");
0069   eccentricity = dqm.book1dHisto("eccentricity", "eccentricity", 200, 0, 1.0, "Eccentricity", "Number of Events");
0070   sigma_inel_NN = dqm.book1dHisto("sigma_inel_NN",
0071                                   "sigma_inel_NN",
0072                                   200,
0073                                   0,
0074                                   10.0,
0075                                   "#sigma{nucleon-nucleon inelastic cross-section}",
0076                                   "Number of Events");
0077 
0078   return;
0079 }
0080 
0081 void BasicHepMCHeavyIonValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0082   ///counters to zero for every event
0083 
0084   ///Gathering the HepMCProduct information
0085   edm::Handle<HepMCProduct> evt;
0086   iEvent.getByToken(hepmcCollectionToken_, evt);
0087 
0088   //Get EVENT
0089   //HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
0090 
0091   const HepMC::HeavyIon* ion = evt->GetEvent()->heavy_ion();
0092 
0093   if (!ion) {
0094     if (QWdebug_)
0095       std::cout << "!!QW!! HeavyIon == null" << std::endl;
0096     return;
0097   }
0098 
0099   double weight = wmanager_.weight(iEvent);
0100   nEvt->Fill(0.5, weight);
0101 
0102   Ncoll_hard->Fill(ion->Ncoll_hard(), weight);
0103   Npart_proj->Fill(ion->Npart_proj(), weight);
0104   Npart_targ->Fill(ion->Npart_targ(), weight);
0105   Ncoll->Fill(ion->Ncoll(), weight);
0106   N_Nwounded_collisions->Fill(ion->N_Nwounded_collisions(), weight);
0107   Nwounded_N_collisions->Fill(ion->Nwounded_N_collisions(), weight);
0108   Nwounded_Nwounded_collisions->Fill(ion->Nwounded_Nwounded_collisions(), weight);
0109   spectator_neutrons->Fill(ion->spectator_neutrons(), weight);
0110   spectator_protons->Fill(ion->spectator_protons(), weight);
0111   impact_parameter->Fill(ion->impact_parameter(), weight);
0112   event_plane_angle->Fill(ion->event_plane_angle(), weight);
0113   eccentricity->Fill(ion->eccentricity(), weight);
0114   sigma_inel_NN->Fill(ion->sigma_inel_NN(), weight);
0115 
0116   //delete myGenEvent;
0117 }  //analyze