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