HiCentralityBiasFilter

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132
// -*- C++ -*-
//
// Package:    HiCentralityBiasFilter
// Class:      HiCentralityBiasFilter
//
/**\class HiCentralityBiasFilter HiCentralityBiasFilter.cc yetkin/HiCentralityBiasFilter/src/HiCentralityBiasFilter.cc

 Description: <one line class summary>

 Implementation:
     <Notes on implementation>
*/
//
// Original Author:  Yetkin Yilmaz
//         Created:  Tue Aug 11 12:42:25 EDT 2009
//
//

// system include files
#include <memory>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDFilter.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "FWCore/Utilities/interface/EDGetToken.h"
#include "FWCore/Utilities/interface/InputTag.h"

#include "FWCore/AbstractServices/interface/RandomNumberGenerator.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"

#include "CLHEP/Random/RandomEngine.h"
#include "HepMC/GenEvent.h"
#include "HepMC/HeavyIon.h"
#include "TF1.h"

using namespace std;

//
// class declaration
//

class HiCentralityBiasFilter : public edm::one::EDFilter<> {
public:
  explicit HiCentralityBiasFilter(const edm::ParameterSet&);
  ~HiCentralityBiasFilter() override = default;

private:
  void beginJob() override;
  bool filter(edm::Event&, const edm::EventSetup&) override;
  void endJob() override;

  const edm::EDGetTokenT<edm::HepMCProduct> hepmcSrc_;
  const std::string func_;
  const std::vector<double> par_;

  edm::Service<edm::RandomNumberGenerator> rng_;

  TF1* fBias_;
};

//
// constants, enums and typedefs
//

//
// static data member definitions
//

//
// constructors and destructor
//
HiCentralityBiasFilter::HiCentralityBiasFilter(const edm::ParameterSet& iConfig)
    : hepmcSrc_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("generatorSmeared"))),
      func_(iConfig.getParameter<string>("function")),
      par_(iConfig.getParameter<vector<double> >("parameters")) {
  // Do what ever initialization is needed
}

//
// member functions
//

// ------------ method called on each new Event  ------------
bool HiCentralityBiasFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
  using namespace edm;

  CLHEP::HepRandomEngine& engine = rng_->getEngine(iEvent.streamID());

  const edm::Handle<HepMCProduct>& mc = iEvent.getHandle(hepmcSrc_);
  const HepMC::GenEvent* evt = mc->GetEvent();

  const HepMC::HeavyIon* hi = evt->heavy_ion();
  if (!hi)
    return false;

  double b = hi->impact_parameter();
  double bound = fBias_->Eval(b);
  double rand = engine.flat();
  if (rand > bound)
    return false;

  return true;
}

// ------------ method called once each job just before starting event loop  ------------
void HiCentralityBiasFilter::beginJob() {
  fBias_ = new TF1("fBias", func_.data(), 0, 20);

  for (size_t ip = 0; ip < par_.size(); ++ip) {
    fBias_->SetParameter(ip, par_[ip]);
  }

  double maxpoint = fBias_->GetMaximum(-0.1, 20);
  if (maxpoint < 0.9)
    throw cms::Exception("HeavyIonCentralityBias")
        << "Input bias function is not optimized. Peak value is " << maxpoint
        << " which is required to be close to 1. Please fix the parameters before production." << endl;
}

// ------------ method called once each job just after ending the event loop  ------------
void HiCentralityBiasFilter::endJob() {}

//define this as a plug-in
DEFINE_FWK_MODULE(HiCentralityBiasFilter);