CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PileUpProducer.cc
Go to the documentation of this file.
5 
7 
9 
16 
17 #include "HepMC/GenEvent.h"
18 #include "TFile.h"
19 #include "TTree.h"
20 #include "TROOT.h"
21 
22 #include <iostream>
23 #include <memory>
24 #include <sys/stat.h>
25 #include <cmath>
26 
28 {
29 
30  // This producer produces an object PileupMixingContent, needed by PileupSummaryInfo
31  produces<PileupMixingContent>();
32  // This producer produces a HepMCProduct, with all pileup vertices/particles
33  produces<edm::HepMCProduct>("PileUpEvents");
34 
35  // Initialize the random number generator service
37  if ( ! rng.isAvailable() ) {
38  throw cms::Exception("Configuration")
39  << "PileUpProducer requires the RandomGeneratorService\n"
40  "which is not present in the configuration file.\n"
41  "You must add the service in the configuration file\n"
42  "or remove the module that requires it";
43  }
44  random = new RandomEngine(&(*rng));
45  gRandom->SetSeed(rng->mySeed());
46 
47  // The pile-up event generation condition
48  const edm::ParameterSet& pu = p.getParameter<edm::ParameterSet>("PileUpSimulator");
49  usePoisson_ = pu.getParameter<bool>("usePoisson");
50  averageNumber_ = 0.;
51  if (usePoisson_) {
52  std::cout << " FastSimulation/PileUpProducer -> poissonian distribution" << std::endl;
53  averageNumber_ = pu.getParameter<double>("averageNumber");
54  } else {//take distribution from configuration
55  dataProbFunctionVar = pu.getParameter<std::vector<int> >("probFunctionVariable");
56  dataProb = pu.getParameter<std::vector<double> >("probValue");
57  varSize = (int) dataProbFunctionVar.size();
58  probSize = (int) dataProb.size();
59  // std::cout << " FastSimulation/PileUpProducer -> varSize = " << varSize << std::endl;
60  // std::cout << " FastSimulation/PileUpProducer -> probSize = " << probSize << std::endl;
61 
62  std::cout << " FastSimulation/PileUpProducer -> distribution from configuration file " << std::endl;
63  if (probSize < varSize){
64  for (int i=0; i<(varSize - probSize); i++) dataProb.push_back(0);
65  edm::LogWarning("") << " The probability function data will be completed with "
66  << (varSize - probSize) <<" values `0';"
67  << " the number of the P(x) data set after adding those 0's is " << dataProb.size();
68  probSize = dataProb.size();
69  }
70  // Create an histogram with the data from the probability function provided by the user
71  int xmin = (int) dataProbFunctionVar[0];
72  int xmax = (int) dataProbFunctionVar[varSize-1]+1; // need upper edge to be one beyond last value
73  int numBins = varSize;
74  std::cout << " FastSimulation/PileUpProducer -> An histogram will be created with " << numBins
75  << " bins in the range ("<< xmin << "," << xmax << ")." << std::endl;
76  hprob = new TH1F("h","Histo from the user's probability function",numBins,xmin,xmax);
77  LogDebug("") << " FastSimulation/PileUpProducer -> Filling histogram with the following data:";
78  for (int j=0; j < numBins ; j++){
79  LogDebug("") << " x = " << dataProbFunctionVar[j ]<< " P(x) = " << dataProb[j];
80  hprob->Fill(dataProbFunctionVar[j]+0.5,dataProb[j]); // assuming integer values for the bins, fill bin centers, not edges
81  }
82 
83  // Check if the histogram is normalized
84  if (((hprob->Integral() - 1) > 1.0e-02) && ((hprob->Integral() - 1) < -1.0e-02)) throw cms::Exception("BadHistoDistribution") << "The histogram should be normalized!" << std::endl;
85 
86  // Get the averageNumber from the histo
87  averageNumber_ = hprob->GetMean();
88  }
89  theFileNames = pu.getUntrackedParameter<std::vector<std::string> >("fileNames");
90  inputFile = pu.getUntrackedParameter<std::string>("inputFile");
92  theFiles.resize(theNumberOfFiles);
93  theTrees.resize(theNumberOfFiles);
94  theBranches.resize(theNumberOfFiles);
95  thePUEvents.resize(theNumberOfFiles);
96  theCurrentEntry.resize(theNumberOfFiles);
97  theCurrentMinBiasEvt.resize(theNumberOfFiles);
98  theNumberOfEntries.resize(theNumberOfFiles);
99  theNumberOfMinBiasEvts.resize(theNumberOfFiles);
100 
101  // Initialize the primary vertex generator
102  const edm::ParameterSet& vtx = p.getParameter<edm::ParameterSet>("VertexGenerator");
103  std::string vtxType = vtx.getParameter<std::string>("type");
104  if ( vtxType == "Gaussian" )
106  else if ( vtxType == "Flat" )
108  else if ( vtxType == "BetaFunc" )
110  else
112 
113 
114 
115  if (averageNumber_ > 0.)
116  {
117  std::cout << " FastSimulation/PileUpProducer -> MinBias events taken from " << theFileNames[0] << " et al., " ;
118  std::cout << " with an average number of events of " << averageNumber_ << std::endl;
119  }
120  else std::cout << " FastSimulation/PileUpProducer -> No pileup " << std::endl;
121 
122 
123 }
124 
126 
127  delete theVertexGenerator;
128  if (hprob) delete hprob;
129 
130 }
131 
133 {
134 
135  gROOT->cd();
136 
137  std::string fullPath;
138 
139  // Read the information from a previous run (to keep reproducibility)
140  bool input = this->read(inputFile);
141  if ( input )
142  std::cout << "***WARNING*** You are reading pile-up information from the file "
143  << inputFile << " created in an earlier run."
144  << std::endl;
145 
146  // Open the file for saving the information of the current run
147  myOutputFile.open ("PileUpOutputFile.txt");
148  myOutputBuffer = 0;
149 
150  // Open the root files
151  for ( unsigned file=0; file<theNumberOfFiles; ++file ) {
152 
153  if (theFileNames[file].find("file:")==0) {
154  fullPath = theFileNames[file].erase(0,5);
155  }
156  else {
157  edm::FileInPath myDataFile("FastSimulation/PileUpProducer/data/"+theFileNames[file]);
158  fullPath = myDataFile.fullPath();
159  }
160  // theFiles[file] = TFile::Open(theFileNames[file].c_str());
161  theFiles[file] = TFile::Open(fullPath.c_str());
162  if ( !theFiles[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
163  << "File " << theFileNames[file] << " " << fullPath << " not found ";
164  //
165  theTrees[file] = (TTree*) theFiles[file]->Get("MinBiasEvents");
166  if ( !theTrees[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
167  << "Tree with name MinBiasEvents not found in " << theFileNames[file];
168  //
169  theBranches[file] = theTrees[file]->GetBranch("puEvent");
170  if ( !theBranches[file] ) throw cms::Exception("FastSimulation/PileUpProducer")
171  << "Branch with name puEvent not found in " << theFileNames[file];
172  //
173  thePUEvents[file] = new PUEvent();
174  theBranches[file]->SetAddress(&thePUEvents[file]);
175  //
176  theNumberOfEntries[file] = theTrees[file]->GetEntries();
177  // Add some randomness (if there was no input file)
178  if ( !input )
180  = (unsigned) (theNumberOfEntries[file] * random->flatShoot());
181 
182  theTrees[file]->GetEntry(theCurrentEntry[file]);
183  unsigned NMinBias = thePUEvents[file]->nMinBias();
184  theNumberOfMinBiasEvts[file] = NMinBias;
185  // Add some randomness (if there was no input file)
186  if ( !input )
188  (unsigned) (theNumberOfMinBiasEvts[file] * random->flatShoot());
189 
190  /*
191  std::cout << "File " << theFileNames[file]
192  << " is opened with " << theNumberOfEntries[file]
193  << " entries and will be read from Entry/Event "
194  << theCurrentEntry[file] << "/" << theCurrentMinBiasEvt[file]
195  << std::endl;
196  */
197  }
198 
199  // Return Loot in the same state as it was when entering.
200  gROOT->cd();
201 
202 }
203 
205 {
206  // Close all local files
207  // Among other things, this allows the TROOT destructor to end up
208  // without crashing, while trying to close these files from outside
209  for ( unsigned file=0; file<theFiles.size(); ++file ) {
210 
211  // std::cout << "Closing " << theFileNames[file] << std::endl;
212  theFiles[file]->Close();
213 
214  }
215 
216  // Close the output file
217  myOutputFile.close();
218 
219  // And return Loot in the same state as it was when entering.
220  gROOT->cd();
221 
222 }
223 
225 {
226 
227  // Create the GenEvent and the HepMCProduct
228  std::auto_ptr<edm::HepMCProduct> pu_product(new edm::HepMCProduct());
229  HepMC::GenEvent* evt = new HepMC::GenEvent();
230 
231  // How many pile-up events?
232  int PUevts; float truePUevts;
233  if (usePoisson_) {
234  PUevts = (int) random->poissonShoot(averageNumber_);
235  truePUevts = (float) averageNumber_;
236  }
237  else {
238  float d = (float) hprob->GetRandom();
239  PUevts = (int) d;
240  truePUevts = d;
241  }
242  // std::cout << "PUevts = " << PUevts << std::endl;
243 
244  // Save this information in the PileupMixingContent object
245  // IMPORTANT: the bunch crossing number is always 0 because FastSim has no out-of-time PU
246  std::auto_ptr< PileupMixingContent > PileupMixing_;
247 
248  std::vector<int> bunchCrossingList;
249  bunchCrossingList.push_back(0);
250 
251  std::vector<int> numInteractionList;
252  numInteractionList.push_back(PUevts);
253 
254  std::vector<float> trueInteractionList;
255  trueInteractionList.push_back(truePUevts);
256 
257  PileupMixing_ = std::auto_ptr< PileupMixingContent >(new PileupMixingContent(bunchCrossingList,numInteractionList,trueInteractionList));
258  iEvent.put(PileupMixing_);
259 
260  // Get N events from random files
261  for ( int ievt=0; ievt<PUevts; ++ievt ) {
262 
263  // Draw a file in a ramdom manner
264  unsigned file = (unsigned) (theNumberOfFiles * random->flatShoot());
265  /*
266  if ( debug )
267  std::cout << "The file chosen for event " << ievt
268  << " is the file number " << file << std::endl;
269  */
270 
271  // Smear the primary vertex and express it in mm (stupid GenEvent convention...)
273  HepMC::FourVector smearedVertex =
274  HepMC::FourVector(theVertexGenerator->X()*10.,
275  theVertexGenerator->Y()*10.,
276  theVertexGenerator->Z()*10.,
277  0.);
278  HepMC::GenVertex* aVertex = new HepMC::GenVertex(smearedVertex);
279  evt->add_vertex(aVertex);
280 
281  // Some rotation around the z axis, for more randomness
282  double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
283  double cAngle = std::cos(theAngle);
284  double sAngle = std::sin(theAngle);
285 
286  /*
287  if ( debug )
288  std::cout << "File chosen : " << file
289  << " Current entry in this file " << theCurrentEntry[file]
290  << " Current minbias in this chunk= " << theCurrentMinBiasEvt[file]
291  << " Total number of minbias in this chunk = " << theNumberOfMinBiasEvts[file] << std::endl;
292  */
293 
294  // theFiles[file]->cd();
295  // gDirectory->ls();
296  // Check we are not either at the end of a minbias bunch
297  // or at the end of a file
298  if ( theCurrentMinBiasEvt[file] == theNumberOfMinBiasEvts[file] ) {
299  // if ( debug ) std::cout << "End of MinBias bunch ! ";
301  // if ( debug) std::cout << "Read the next entry " << theCurrentEntry[file] << std::endl;
303  if ( theCurrentEntry[file] == theNumberOfEntries[file] ) {
304  theCurrentEntry[file] = 0;
305  // if ( debug ) std::cout << "End of file - Rewind! " << std::endl;
306  }
307  // if ( debug ) std::cout << "The PUEvent is reset ... ";
308  thePUEvents[file]->reset();
309  unsigned myEntry = theCurrentEntry[file];
310  /*
311  if ( debug ) std::cout << "The new entry " << myEntry
312  << " is read ... in TTree " << theTrees[file] << " ";
313  */
314  theTrees[file]->GetEntry(myEntry);
315  /*
316  if ( debug )
317  std::cout << "The number of interactions in the new entry is ... ";
318  */
320  // if ( debug ) std::cout << theNumberOfMinBiasEvts[file] << std::endl;
321  }
322 
323  // Read a minbias event chunk
324  const PUEvent::PUMinBiasEvt& aMinBiasEvt
325  = thePUEvents[file]->thePUMinBiasEvts()[theCurrentMinBiasEvt[file]];
326 
327  // Find corresponding particles
328  unsigned firstTrack = aMinBiasEvt.first;
329  unsigned trackSize = firstTrack + aMinBiasEvt.size;
330  /*
331  if ( debug ) std::cout << "First and last+1 tracks are "
332  << firstTrack << " " << trackSize << std::endl;
333  */
334 
335  // Loop on particles
336  for ( unsigned iTrack=firstTrack; iTrack<trackSize; ++iTrack ) {
337 
338  const PUEvent::PUParticle& aParticle
339  = thePUEvents[file]->thePUParticles()[iTrack];
340  /*
341  if ( debug)
342  std::cout << "Track " << iTrack
343  << " id/px/py/pz/mass "
344  << aParticle.id << " "
345  << aParticle.px << " "
346  << aParticle.py << " "
347  << aParticle.pz << " "
348  << aParticle.mass << " " << std::endl;
349  */
350 
351  // Create a FourVector, with rotation
352  double energy = std::sqrt( aParticle.px*aParticle.px
353  + aParticle.py*aParticle.py
354  + aParticle.pz*aParticle.pz
355  + aParticle.mass*aParticle.mass );
356 
357  HepMC::FourVector myPart(cAngle * aParticle.px + sAngle * aParticle.py,
358  -sAngle * aParticle.px + cAngle * aParticle.py,
359  aParticle.pz, energy);
360 
361  // Add a GenParticle
362  HepMC::GenParticle* aGenParticle = new HepMC::GenParticle(myPart,aParticle.id);
363  aVertex->add_particle_out(aGenParticle);
364 
365  }
366  // End of particle loop
367 
368  // Increment for next time
370 
371  }
372  // End of pile-up event loop
373 
374  // evt->print();
375 
376  // Fill the HepMCProduct from the GenEvent
377  if ( evt ) {
378  pu_product->addHepMCData( evt );
379  // Boost in case of beam crossing angle
380  TMatrixD* boost = theVertexGenerator->boost();
381  if ( boost ) pu_product->boostToLab(boost,"momentum");
382  }
383 
384  // Put the HepMCProduct onto the event
385  iEvent.put(pu_product,"PileUpEvents");
386  // delete evt;
387 
388  // Save the current location in each pile-up event files
389  this->save();
390 
391 }
392 
393 void
395 
396  // Size of buffer
397  ++myOutputBuffer;
398 
399  // Periodically close the current file and open a new one
400  if ( myOutputBuffer/1000*1000 == myOutputBuffer ) {
401  myOutputFile.close();
402  myOutputFile.open ("PileUpOutputFile.txt");
403  // myOutputFile.seekp(0); // No need to rewind in that case
404  }
405 
406  // Save the current position to file
407  myOutputFile.write((const char*)(&theCurrentEntry.front()),
408  theCurrentEntry.size()*sizeof(unsigned));
409  myOutputFile.write((const char*)&theCurrentMinBiasEvt.front(),
410  theCurrentMinBiasEvt.size()*sizeof(unsigned));
411  myOutputFile.flush();
412 
413 }
414 
415 bool
417 
418  ifstream myInputFile;
419  struct stat results;
420  unsigned size1 = theCurrentEntry.size()*sizeof(unsigned);
421  unsigned size2 = theCurrentMinBiasEvt.size()*sizeof(unsigned);
422  unsigned size = 0;
423 
424 
425  // Open the file (if any)
426  myInputFile.open (inputFile.c_str());
427  if ( myInputFile.is_open() ) {
428 
429  // Get the size of the file
430  if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
431  else return false; // Something is wrong with that file !
432 
433  // Position the pointer just before the last record
434  myInputFile.seekg(size-size1-size2);
435 
436  // Read the information
437  myInputFile.read((char*)(&theCurrentEntry.front()),size1);
438  myInputFile.read((char*)&theCurrentMinBiasEvt.front(),size2);
439  myInputFile.close();
440 
441  return true;
442 
443  }
444 
445  return false;
446 
447 }
448 
#define LogDebug(id)
PileUpProducer(edm::ParameterSet const &p)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< unsigned > theNumberOfMinBiasEvts
std::vector< std::string > theFileNames
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool read(std::string inputFile)
Read former minbias configuration (from previous run)
std::vector< int > dataProbFunctionVar
virtual void endRun()
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< unsigned > theNumberOfEntries
PrimaryVertexGenerator * theVertexGenerator
std::ofstream myOutputFile
int iEvent
Definition: GenABIO.cc:243
virtual void produce(edm::Event &e, const edm::EventSetup &c)
virtual void beginRun(edm::Run &, edm::EventSetup const &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
T sqrt(T t)
Definition: SSEVec.h:46
unsigned int poissonShoot(double mean) const
Definition: RandomEngine.h:44
const RandomEngine * random
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< unsigned > theCurrentMinBiasEvt
int j
Definition: DBlmapReader.cc:9
std::vector< PUEvent * > thePUEvents
std::string inputFile
virtual void generate()=0
Generation process (to be implemented)
std::vector< double > dataProb
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
std::vector< unsigned > theCurrentEntry
double averageNumber_
std::vector< TTree * > theTrees
virtual ~PileUpProducer()
Definition: PUEvent.h:6
tuple cout
Definition: gather_cfg.py:121
std::string fullPath() const
Definition: FileInPath.cc:171
unsigned theNumberOfFiles
unsigned myOutputBuffer
tuple size
Write out results.
Definition: Run.h:33
std::vector< TFile * > theFiles
void save()
Save current minbias configuration (for later use)
std::vector< TBranch * > theBranches