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/HepMCProduct/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)
00028 {
00029
00030
00031 produces<edm::HepMCProduct>("PileUpEvents");
00032
00033
00034 edm::Service<edm::RandomNumberGenerator> rng;
00035 if ( ! rng.isAvailable() ) {
00036 throw cms::Exception("Configuration")
00037 << "PileUpProducer requires the RandomGeneratorService\n"
00038 "which is not present in the configuration file.\n"
00039 "You must add the service in the configuration file\n"
00040 "or remove the module that requires it";
00041 }
00042 random = new RandomEngine(&(*rng));
00043
00044
00045 const edm::ParameterSet& pu = p.getParameter<edm::ParameterSet>("PileUpSimulator");
00046 averageNumber_ = pu.getParameter<double>("averageNumber");
00047 theFileNames = pu.getUntrackedParameter<std::vector<std::string> >("fileNames");
00048 inputFile = pu.getUntrackedParameter<std::string>("inputFile");
00049 theNumberOfFiles = theFileNames.size();
00050 theFiles.resize(theNumberOfFiles);
00051 theTrees.resize(theNumberOfFiles);
00052 theBranches.resize(theNumberOfFiles);
00053 thePUEvents.resize(theNumberOfFiles);
00054 theCurrentEntry.resize(theNumberOfFiles);
00055 theCurrentMinBiasEvt.resize(theNumberOfFiles);
00056 theNumberOfEntries.resize(theNumberOfFiles);
00057 theNumberOfMinBiasEvts.resize(theNumberOfFiles);
00058
00059
00060 const edm::ParameterSet& vtx = p.getParameter<edm::ParameterSet>("VertexGenerator");
00061 std::string vtxType = vtx.getParameter<std::string>("type");
00062 if ( vtxType == "Gaussian" )
00063 theVertexGenerator = new GaussianPrimaryVertexGenerator(vtx,random);
00064 else if ( vtxType == "Flat" )
00065 theVertexGenerator = new FlatPrimaryVertexGenerator(vtx,random);
00066 else if ( vtxType == "BetaFunc" )
00067 theVertexGenerator = new BetaFuncPrimaryVertexGenerator(vtx,random);
00068 else
00069 theVertexGenerator = new NoPrimaryVertexGenerator();
00070
00071 }
00072
00073 PileUpProducer::~PileUpProducer() {
00074
00075 delete theVertexGenerator;
00076
00077 }
00078
00079 void PileUpProducer::beginJob(const edm::EventSetup & es)
00080 {
00081
00082 std::cout << " PileUpProducer initializing " << std::endl;
00083 gROOT->cd();
00084
00085 std::string fullPath;
00086
00087
00088 bool input = this->read(inputFile);
00089 if ( input )
00090 std::cout << "***WARNING*** You are reading pile-up information from the file "
00091 << inputFile << " created in an earlier run."
00092 << std::endl;
00093
00094
00095 myOutputFile.open ("PileUpOutputFile.txt");
00096 myOutputBuffer = 0;
00097
00098
00099 std::cout << "Opening minimum-bias event files ... " << std::endl;
00100 for ( unsigned file=0; file<theNumberOfFiles; ++file ) {
00101
00102 edm::FileInPath myDataFile("FastSimulation/PileUpProducer/data/"+theFileNames[file]);
00103 fullPath = myDataFile.fullPath();
00104
00105 theFiles[file] = TFile::Open(fullPath.c_str());
00106 if ( !theFiles[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00107 << "File " << theFileNames[file] << " " << fullPath << " not found ";
00108
00109 theTrees[file] = (TTree*) theFiles[file]->Get("MinBiasEvents");
00110 if ( !theTrees[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00111 << "Tree with name MinBiasEvents not found in " << theFileNames[file];
00112
00113 theBranches[file] = theTrees[file]->GetBranch("puEvent");
00114 if ( !theBranches[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
00115 << "Branch with name puEvent not found in " << theFileNames[file];
00116
00117 thePUEvents[file] = new PUEvent();
00118 theBranches[file]->SetAddress(&thePUEvents[file]);
00119
00120 theNumberOfEntries[file] = theTrees[file]->GetEntries();
00121
00122 if ( !input )
00123 theCurrentEntry[file]
00124 = (unsigned) (theNumberOfEntries[file] * random->flatShoot());
00125
00126 theTrees[file]->GetEntry(theCurrentEntry[file]);
00127 unsigned NMinBias = thePUEvents[file]->nMinBias();
00128 theNumberOfMinBiasEvts[file] = NMinBias;
00129
00130 if ( !input )
00131 theCurrentMinBiasEvt[file] =
00132 (unsigned) (theNumberOfMinBiasEvts[file] * random->flatShoot());
00133
00134
00135
00136
00137
00138
00139
00140
00141 }
00142
00143
00144 gROOT->cd();
00145
00146 }
00147
00148 void PileUpProducer::endJob()
00149 {
00150 std::cout << " PileUpProducer terminating " << std::endl;
00151
00152
00153
00154 std::cout << "Closing minimum-bias event files... " << std::endl;
00155 for ( unsigned file=0; file<theFiles.size(); ++file ) {
00156
00157
00158 theFiles[file]->Close();
00159
00160 }
00161
00162
00163 myOutputFile.close();
00164
00165
00166 gROOT->cd();
00167
00168 }
00169
00170 void PileUpProducer::produce(edm::Event & iEvent, const edm::EventSetup & es)
00171 {
00172
00173
00174 std::auto_ptr<edm::HepMCProduct> pu_product(new edm::HepMCProduct());
00175 HepMC::GenEvent* evt = new HepMC::GenEvent();
00176
00177
00178 int PUevts = (int) random->poissonShoot(averageNumber_);
00179
00180
00181 for ( int ievt=0; ievt<PUevts; ++ievt ) {
00182
00183
00184 unsigned file = (unsigned) (theNumberOfFiles * random->flatShoot());
00185
00186
00187
00188
00189
00190
00191
00192 theVertexGenerator->generate();
00193 HepMC::FourVector smearedVertex =
00194 HepMC::FourVector(theVertexGenerator->X()*10.,
00195 theVertexGenerator->Y()*10.,
00196 theVertexGenerator->Z()*10.,
00197 0.);
00198 HepMC::GenVertex* aVertex = new HepMC::GenVertex(smearedVertex);
00199 evt->add_vertex(aVertex);
00200
00201
00202 double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
00203 double cAngle = std::cos(theAngle);
00204 double sAngle = std::sin(theAngle);
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218 if ( theCurrentMinBiasEvt[file] == theNumberOfMinBiasEvts[file] ) {
00219
00220 ++theCurrentEntry[file];
00221
00222 theCurrentMinBiasEvt[file] = 0;
00223 if ( theCurrentEntry[file] == theNumberOfEntries[file] ) {
00224 theCurrentEntry[file] = 0;
00225
00226 }
00227
00228 thePUEvents[file]->reset();
00229 unsigned myEntry = theCurrentEntry[file];
00230
00231
00232
00233
00234 theTrees[file]->GetEntry(myEntry);
00235
00236
00237
00238
00239 theNumberOfMinBiasEvts[file] = thePUEvents[file]->nMinBias();
00240
00241 }
00242
00243
00244 const PUEvent::PUMinBiasEvt& aMinBiasEvt
00245 = thePUEvents[file]->thePUMinBiasEvts()[theCurrentMinBiasEvt[file]];
00246
00247
00248 unsigned firstTrack = aMinBiasEvt.first;
00249 unsigned trackSize = firstTrack + aMinBiasEvt.size;
00250
00251
00252
00253
00254
00255
00256 for ( unsigned iTrack=firstTrack; iTrack<trackSize; ++iTrack ) {
00257
00258 const PUEvent::PUParticle& aParticle
00259 = thePUEvents[file]->thePUParticles()[iTrack];
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272 double energy = std::sqrt( aParticle.px*aParticle.px
00273 + aParticle.py*aParticle.py
00274 + aParticle.pz*aParticle.pz
00275 + aParticle.mass*aParticle.mass );
00276
00277 HepMC::FourVector myPart(cAngle * aParticle.px + sAngle * aParticle.py,
00278 -sAngle * aParticle.px + cAngle * aParticle.py,
00279 aParticle.pz, energy);
00280
00281
00282 HepMC::GenParticle* aGenParticle = new HepMC::GenParticle(myPart,aParticle.id);
00283 aVertex->add_particle_out(aGenParticle);
00284
00285 }
00286
00287
00288
00289 ++theCurrentMinBiasEvt[file];
00290
00291 }
00292
00293
00294
00295
00296
00297 if ( evt ) {
00298 pu_product->addHepMCData( evt );
00299
00300 TMatrixD* boost = theVertexGenerator->boost();
00301 if ( boost ) pu_product->boostToLab(boost,"momentum");
00302 }
00303
00304
00305 iEvent.put(pu_product,"PileUpEvents");
00306
00307
00308
00309 this->save();
00310
00311 }
00312
00313 void
00314 PileUpProducer::save() {
00315
00316
00317 ++myOutputBuffer;
00318
00319
00320 if ( myOutputBuffer/1000*1000 == myOutputBuffer ) {
00321 myOutputFile.close();
00322 myOutputFile.open ("PileUpOutputFile.txt");
00323
00324 }
00325
00326
00327 myOutputFile.write((const char*)(&theCurrentEntry.front()),
00328 theCurrentEntry.size()*sizeof(unsigned));
00329 myOutputFile.write((const char*)&theCurrentMinBiasEvt.front(),
00330 theCurrentMinBiasEvt.size()*sizeof(unsigned));
00331 myOutputFile.flush();
00332
00333 }
00334
00335 bool
00336 PileUpProducer::read(std::string inputFile) {
00337
00338 ifstream myInputFile;
00339 struct stat results;
00340 unsigned size1 = theCurrentEntry.size()*sizeof(unsigned);
00341 unsigned size2 = theCurrentMinBiasEvt.size()*sizeof(unsigned);
00342 unsigned size = 0;
00343
00344
00345
00346 myInputFile.open (inputFile.c_str());
00347 if ( myInputFile.is_open() ) {
00348
00349
00350 if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
00351 else return false;
00352
00353
00354 myInputFile.seekg(size-size1-size2);
00355
00356
00357 myInputFile.read((char*)(&theCurrentEntry.front()),size1);
00358 myInputFile.read((char*)&theCurrentMinBiasEvt.front(),size2);
00359 myInputFile.close();
00360
00361 return true;
00362
00363 }
00364
00365 return false;
00366
00367 }
00368
00369 DEFINE_FWK_MODULE(PileUpProducer);