Go to the documentation of this file.00001 #include "FWCore/Framework/interface/Event.h"
00002 #include "FWCore/Framework/interface/MakerMacros.h"
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00005
00006 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00007
00008 #include "FastSimDataFormats/PileUpEvents/interface/PUEvent.h"
00009
00010 #include "FastSimulation/PileUpProducer/plugins/PileUpProducer.h"
00011 #include "FastSimulation/Event/interface/BetaFuncPrimaryVertexGenerator.h"
00012 #include "FastSimulation/Event/interface/GaussianPrimaryVertexGenerator.h"
00013 #include "FastSimulation/Event/interface/FlatPrimaryVertexGenerator.h"
00014 #include "FastSimulation/Event/interface/NoPrimaryVertexGenerator.h"
00015 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00016
00017 #include "HepMC/GenEvent.h"
00018 #include "TFile.h"
00019 #include "TTree.h"
00020 #include "TROOT.h"
00021
00022 #include <iostream>
00023 #include <memory>
00024 #include <sys/stat.h>
00025 #include <cmath>
00026
00027 PileUpProducer::PileUpProducer(edm::ParameterSet const & p) : hprob(0)
00028 {
00029
00030
00031 produces<PileupMixingContent>();
00032
00033 produces<edm::HepMCProduct>("PileUpEvents");
00034
00035
00036 edm::Service<edm::RandomNumberGenerator> rng;
00037 if ( ! rng.isAvailable() ) {
00038 throw cms::Exception("Configuration")
00039 << "PileUpProducer requires the RandomGeneratorService\n"
00040 "which is not present in the configuration file.\n"
00041 "You must add the service in the configuration file\n"
00042 "or remove the module that requires it";
00043 }
00044 random = new RandomEngine(&(*rng));
00045
00046
00047 const edm::ParameterSet& pu = p.getParameter<edm::ParameterSet>("PileUpSimulator");
00048 usePoisson_ = pu.getParameter<bool>("usePoisson");
00049 averageNumber_ = 0.;
00050 if (usePoisson_) {
00051 std::cout << " FastSimulation/PileUpProducer -> poissonian distribution" << std::endl;
00052 averageNumber_ = pu.getParameter<double>("averageNumber");
00053 } else {
00054 dataProbFunctionVar = pu.getParameter<std::vector<int> >("probFunctionVariable");
00055 dataProb = pu.getParameter<std::vector<double> >("probValue");
00056 varSize = (int) dataProbFunctionVar.size();
00057 probSize = (int) dataProb.size();
00058
00059
00060
00061 std::cout << " FastSimulation/PileUpProducer -> distribution from configuration file " << std::endl;
00062 if (probSize < varSize){
00063 for (int i=0; i<(varSize - probSize); i++) dataProb.push_back(0);
00064 edm::LogWarning("") << " The probability function data will be completed with "
00065 << (varSize - probSize) <<" values `0';"
00066 << " the number of the P(x) data set after adding those 0's is " << dataProb.size();
00067 probSize = dataProb.size();
00068 }
00069
00070 int xmin = (int) dataProbFunctionVar[0];
00071 int xmax = (int) dataProbFunctionVar[varSize-1]+1;
00072 int numBins = varSize;
00073 std::cout << " FastSimulation/PileUpProducer -> An histogram will be created with " << numBins
00074 << " bins in the range ("<< xmin << "," << xmax << ")." << std::endl;
00075 hprob = new TH1F("h","Histo from the user's probability function",numBins,xmin,xmax);
00076 LogDebug("") << " FastSimulation/PileUpProducer -> Filling histogram with the following data:";
00077 for (int j=0; j < numBins ; j++){
00078 LogDebug("") << " x = " << dataProbFunctionVar[j ]<< " P(x) = " << dataProb[j];
00079 hprob->Fill(dataProbFunctionVar[j]+0.5,dataProb[j]);
00080 }
00081
00082
00083 if (((hprob->Integral() - 1) > 1.0e-02) && ((hprob->Integral() - 1) < -1.0e-02)) throw cms::Exception("BadHistoDistribution") << "The histogram should be normalized!" << std::endl;
00084
00085
00086 averageNumber_ = hprob->GetMean();
00087 }
00088 theFileNames = pu.getUntrackedParameter<std::vector<std::string> >("fileNames");
00089 inputFile = pu.getUntrackedParameter<std::string>("inputFile");
00090 theNumberOfFiles = theFileNames.size();
00091 theFiles.resize(theNumberOfFiles);
00092 theTrees.resize(theNumberOfFiles);
00093 theBranches.resize(theNumberOfFiles);
00094 thePUEvents.resize(theNumberOfFiles);
00095 theCurrentEntry.resize(theNumberOfFiles);
00096 theCurrentMinBiasEvt.resize(theNumberOfFiles);
00097 theNumberOfEntries.resize(theNumberOfFiles);
00098 theNumberOfMinBiasEvts.resize(theNumberOfFiles);
00099
00100
00101 const edm::ParameterSet& vtx = p.getParameter<edm::ParameterSet>("VertexGenerator");
00102 std::string vtxType = vtx.getParameter<std::string>("type");
00103 if ( vtxType == "Gaussian" )
00104 theVertexGenerator = new GaussianPrimaryVertexGenerator(vtx,random);
00105 else if ( vtxType == "Flat" )
00106 theVertexGenerator = new FlatPrimaryVertexGenerator(vtx,random);
00107 else if ( vtxType == "BetaFunc" )
00108 theVertexGenerator = new BetaFuncPrimaryVertexGenerator(vtx,random);
00109 else
00110 theVertexGenerator = new NoPrimaryVertexGenerator();
00111
00112
00113
00114 if (averageNumber_ > 0.)
00115 {
00116 std::cout << " FastSimulation/PileUpProducer -> MinBias events taken from " << theFileNames[0] << " et al., " ;
00117 std::cout << " with an average number of events of " << averageNumber_ << std::endl;
00118 }
00119 else std::cout << " FastSimulation/PileUpProducer -> No pileup " << std::endl;
00120
00121
00122 }
00123
00124 PileUpProducer::~PileUpProducer() {
00125
00126 delete theVertexGenerator;
00127 if (hprob) delete hprob;
00128
00129 }
00130
00131 void PileUpProducer::beginRun(edm::Run & run, edm::EventSetup const& es)
00132 {
00133
00134 gROOT->cd();
00135
00136 std::string fullPath;
00137
00138
00139 bool input = this->read(inputFile);
00140 if ( input )
00141 std::cout << "***WARNING*** You are reading pile-up information from the file "
00142 << inputFile << " created in an earlier run."
00143 << std::endl;
00144
00145
00146 myOutputFile.open ("PileUpOutputFile.txt");
00147 myOutputBuffer = 0;
00148
00149
00150 for ( unsigned file=0; file<theNumberOfFiles; ++file ) {
00151
00152 if (theFileNames[file].find("file:")==0) {
00153 fullPath = theFileNames[file].erase(0,5);
00154 }
00155 else {
00156 edm::FileInPath myDataFile("FastSimulation/PileUpProducer/data/"+theFileNames[file]);
00157 fullPath = myDataFile.fullPath();
00158 }
00159
00160 theFiles[file] = TFile::Open(fullPath.c_str());
00161 if ( !theFiles[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00162 << "File " << theFileNames[file] << " " << fullPath << " not found ";
00163
00164 theTrees[file] = (TTree*) theFiles[file]->Get("MinBiasEvents");
00165 if ( !theTrees[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00166 << "Tree with name MinBiasEvents not found in " << theFileNames[file];
00167
00168 theBranches[file] = theTrees[file]->GetBranch("puEvent");
00169 if ( !theBranches[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00170 << "Branch with name puEvent not found in " << theFileNames[file];
00171
00172 thePUEvents[file] = new PUEvent();
00173 theBranches[file]->SetAddress(&thePUEvents[file]);
00174
00175 theNumberOfEntries[file] = theTrees[file]->GetEntries();
00176
00177 if ( !input )
00178 theCurrentEntry[file]
00179 = (unsigned) (theNumberOfEntries[file] * random->flatShoot());
00180
00181 theTrees[file]->GetEntry(theCurrentEntry[file]);
00182 unsigned NMinBias = thePUEvents[file]->nMinBias();
00183 theNumberOfMinBiasEvts[file] = NMinBias;
00184
00185 if ( !input )
00186 theCurrentMinBiasEvt[file] =
00187 (unsigned) (theNumberOfMinBiasEvts[file] * random->flatShoot());
00188
00189
00190
00191
00192
00193
00194
00195
00196 }
00197
00198
00199 gROOT->cd();
00200
00201 }
00202
00203 void PileUpProducer::endRun()
00204 {
00205
00206
00207
00208 for ( unsigned file=0; file<theFiles.size(); ++file ) {
00209
00210
00211 theFiles[file]->Close();
00212
00213 }
00214
00215
00216 myOutputFile.close();
00217
00218
00219 gROOT->cd();
00220
00221 }
00222
00223 void PileUpProducer::produce(edm::Event & iEvent, const edm::EventSetup & es)
00224 {
00225
00226
00227 std::auto_ptr<edm::HepMCProduct> pu_product(new edm::HepMCProduct());
00228 HepMC::GenEvent* evt = new HepMC::GenEvent();
00229
00230
00231 int PUevts; float truePUevts;
00232 if (usePoisson_) {
00233 PUevts = (int) random->poissonShoot(averageNumber_);
00234 truePUevts = (float) averageNumber_;
00235 }
00236 else {
00237 float d = (float) hprob->GetRandom();
00238 PUevts = (int) d;
00239 truePUevts = d;
00240 }
00241
00242
00243
00244
00245 std::auto_ptr< PileupMixingContent > PileupMixing_;
00246
00247 std::vector<int> bunchCrossingList;
00248 bunchCrossingList.push_back(0);
00249
00250 std::vector<int> numInteractionList;
00251 numInteractionList.push_back(PUevts);
00252
00253 std::vector<float> trueInteractionList;
00254 trueInteractionList.push_back(truePUevts);
00255
00256 PileupMixing_ = std::auto_ptr< PileupMixingContent >(new PileupMixingContent(bunchCrossingList,numInteractionList,trueInteractionList));
00257 iEvent.put(PileupMixing_);
00258
00259
00260 for ( int ievt=0; ievt<PUevts; ++ievt ) {
00261
00262
00263 unsigned file = (unsigned) (theNumberOfFiles * random->flatShoot());
00264
00265
00266
00267
00268
00269
00270
00271 theVertexGenerator->generate();
00272 HepMC::FourVector smearedVertex =
00273 HepMC::FourVector(theVertexGenerator->X()*10.,
00274 theVertexGenerator->Y()*10.,
00275 theVertexGenerator->Z()*10.,
00276 0.);
00277 HepMC::GenVertex* aVertex = new HepMC::GenVertex(smearedVertex);
00278 evt->add_vertex(aVertex);
00279
00280
00281 double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
00282 double cAngle = std::cos(theAngle);
00283 double sAngle = std::sin(theAngle);
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297 if ( theCurrentMinBiasEvt[file] == theNumberOfMinBiasEvts[file] ) {
00298
00299 ++theCurrentEntry[file];
00300
00301 theCurrentMinBiasEvt[file] = 0;
00302 if ( theCurrentEntry[file] == theNumberOfEntries[file] ) {
00303 theCurrentEntry[file] = 0;
00304
00305 }
00306
00307 thePUEvents[file]->reset();
00308 unsigned myEntry = theCurrentEntry[file];
00309
00310
00311
00312
00313 theTrees[file]->GetEntry(myEntry);
00314
00315
00316
00317
00318 theNumberOfMinBiasEvts[file] = thePUEvents[file]->nMinBias();
00319
00320 }
00321
00322
00323 const PUEvent::PUMinBiasEvt& aMinBiasEvt
00324 = thePUEvents[file]->thePUMinBiasEvts()[theCurrentMinBiasEvt[file]];
00325
00326
00327 unsigned firstTrack = aMinBiasEvt.first;
00328 unsigned trackSize = firstTrack + aMinBiasEvt.size;
00329
00330
00331
00332
00333
00334
00335 for ( unsigned iTrack=firstTrack; iTrack<trackSize; ++iTrack ) {
00336
00337 const PUEvent::PUParticle& aParticle
00338 = thePUEvents[file]->thePUParticles()[iTrack];
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351 double energy = std::sqrt( aParticle.px*aParticle.px
00352 + aParticle.py*aParticle.py
00353 + aParticle.pz*aParticle.pz
00354 + aParticle.mass*aParticle.mass );
00355
00356 HepMC::FourVector myPart(cAngle * aParticle.px + sAngle * aParticle.py,
00357 -sAngle * aParticle.px + cAngle * aParticle.py,
00358 aParticle.pz, energy);
00359
00360
00361 HepMC::GenParticle* aGenParticle = new HepMC::GenParticle(myPart,aParticle.id);
00362 aVertex->add_particle_out(aGenParticle);
00363
00364 }
00365
00366
00367
00368 ++theCurrentMinBiasEvt[file];
00369
00370 }
00371
00372
00373
00374
00375
00376 if ( evt ) {
00377 pu_product->addHepMCData( evt );
00378
00379 TMatrixD* boost = theVertexGenerator->boost();
00380 if ( boost ) pu_product->boostToLab(boost,"momentum");
00381 }
00382
00383
00384 iEvent.put(pu_product,"PileUpEvents");
00385
00386
00387
00388 this->save();
00389
00390 }
00391
00392 void
00393 PileUpProducer::save() {
00394
00395
00396 ++myOutputBuffer;
00397
00398
00399 if ( myOutputBuffer/1000*1000 == myOutputBuffer ) {
00400 myOutputFile.close();
00401 myOutputFile.open ("PileUpOutputFile.txt");
00402
00403 }
00404
00405
00406 myOutputFile.write((const char*)(&theCurrentEntry.front()),
00407 theCurrentEntry.size()*sizeof(unsigned));
00408 myOutputFile.write((const char*)&theCurrentMinBiasEvt.front(),
00409 theCurrentMinBiasEvt.size()*sizeof(unsigned));
00410 myOutputFile.flush();
00411
00412 }
00413
00414 bool
00415 PileUpProducer::read(std::string inputFile) {
00416
00417 ifstream myInputFile;
00418 struct stat results;
00419 unsigned size1 = theCurrentEntry.size()*sizeof(unsigned);
00420 unsigned size2 = theCurrentMinBiasEvt.size()*sizeof(unsigned);
00421 unsigned size = 0;
00422
00423
00424
00425 myInputFile.open (inputFile.c_str());
00426 if ( myInputFile.is_open() ) {
00427
00428
00429 if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
00430 else return false;
00431
00432
00433 myInputFile.seekg(size-size1-size2);
00434
00435
00436 myInputFile.read((char*)(&theCurrentEntry.front()),size1);
00437 myInputFile.read((char*)&theCurrentMinBiasEvt.front(),size2);
00438 myInputFile.close();
00439
00440 return true;
00441
00442 }
00443
00444 return false;
00445
00446 }
00447
00448 DEFINE_FWK_MODULE(PileUpProducer);