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