CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TreeWriterForEcalCorrection.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
7 
10 
12 
15 
17 
20 
21 #include "TTree.h"
22 
23 
25  public:
28 
29  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
30 
31 
32  private:
33  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
34 
36  TTree* tree;
38 };
39 
41 {
42  tree = file->make<TTree>("responseTree", "same info as 3dhisto");
43  tree->Branch( "e", &tree_e, "e/F");
44  tree->Branch( "eta", &tree_eta, "eta/F");
45  tree->Branch( "r", &tree_response, "r/F");
46 }
47 
48 void
50 {
51  // get generated particles
53  iEvent.getByLabel("genParticles", "", GenParticles);
54 
55  // As this module is intended for single particle guns, there should be
56  // exactly one generated particle.
57  if( GenParticles->size() != 1 ) {
58  throw cms::Exception("MismatchedInputFiles")
59  << "Intended for particle guns only\n";
60  }
61 
62  // I assume here that the tracker simulation is disabled and no vertex
63  // smearing is done, so that the generated particles is exactly the same
64  // particle as the particle hitting the entrace of the ECAL.
65  reco::GenParticle gen = GenParticles->at(0);
66  float genE = gen.energy();
67  float genEta = gen.eta();
68  // genEta is positive for my photongun sample per definition.
69  // If not, you could take the absolute value here.
70 
71  // get sim hits
75 
76  // Finds out automatically, if this is fullsim or fastsim
77  bool isFastSim = iEvent.getByLabel( "famosSimHits", "EcalHitsEB", SimHitsEB );
78  if( isFastSim ) {
79  iEvent.getByLabel( "famosSimHits", "EcalHitsEE", SimHitsEE );
80  iEvent.getByLabel( "famosSimHits", "EcalHitsES", SimHitsES );
81  } else {
82  iEvent.getByLabel( "g4SimHits", "EcalHitsEB", SimHitsEB );
83  iEvent.getByLabel( "g4SimHits", "EcalHitsEE", SimHitsEE );
84  iEvent.getByLabel( "g4SimHits", "EcalHitsES", SimHitsES );
85  }
86 
87  // merge them into one single vector
88  auto SimHits = *SimHitsEB;
89  SimHits.insert( SimHits.end(), SimHitsEE->begin(), SimHitsEE->end() );
90  SimHits.insert( SimHits.end(), SimHitsES->begin(), SimHitsES->end() );
91 
92  // As we only had one generated particle (and hopefully no pileup),
93  // the total energy is due to the generated particle only
94  float energyTotal = 0;
95  for( auto const& Hit : SimHits ) {
96  energyTotal += Hit.energy();
97  }
98 
99  tree_e = genE;
100  tree_eta = genEta;
101  tree_response = energyTotal/genE;
102  tree->Fill();
103 
104 }
105 
106 void
109  descriptions.add( "ecalScaleFactorCalculator", desc );
110 }
111 
edm::Service< TFileService > file
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
virtual double eta() const
momentum pseudorapidity
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:230
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TreeWriterForEcalCorrection(const edm::ParameterSet &)