LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
CORSIKAGen_module.cc
Go to the documentation of this file.
1 #ifndef EVGEN_CORSIKAGen_H
8 #define EVGEN_CORSIKAGen_H
9 
10 // ROOT includes
11 #include "TRandom3.h"
12 #include "TDatabasePDG.h"
13 #include "TString.h"
14 #include "TSystem.h" //need BaseName and DirName
15 
16 // Framework includes
20 #include "fhiclcpp/ParameterSet.h"
28 
29 // art extensions
31 
32 // larsoft includes
40 
41 #include <sqlite3.h>
42 #include "CLHEP/Random/RandFlat.h"
43 #include "CLHEP/Random/RandPoissonQ.h"
44 #include "ifdh.h" //to handle flux files
45 
46 namespace evgen {
47 
49  class CORSIKAGen : public art::EDProducer {
50  public:
51  explicit CORSIKAGen(fhicl::ParameterSet const& pset);
52  virtual ~CORSIKAGen();
53 
54  void produce(art::Event& evt);
55  void beginJob();
56  void beginRun(art::Run& run);
57  void reconfigure(fhicl::ParameterSet const& p);
58 
59 
60  private:
61  void openDBs();
62  void populateNShowers();
63  void populateTOffset();
64  void GetSample(simb::MCTruth&);
65  double wrapvar( const double var, const double low, const double high);
66  double wrapvarBoxNo( const double var, const double low, const double high, int& boxno);
67  void ProjectToBoxEdge(const double xyz[],
68  const double dxyz[],
69  const double xlo,
70  const double xhi,
71  const double ylo,
72  const double yhi,
73  const double zlo,
74  const double zhi,
75  double xyzout[]);
76 
77  int fShowerInputs=0;
78  std::vector<double> fNShowersPerEvent;
79  std::vector<int> fMaxShowers; //< Max number of showers to query, one per showerinput
80  double fShowerBounds[6]={0.,0.,0.,0.,0.,0.};
81  double fToffset_corsika=0.;
82  ifdh_ns::ifdh* fIFDH=0;
83 
84  //fcl parameters
85  double fProjectToHeight=0.;
86  std::vector< std::string > fShowerInputFiles;
87  std::vector< double > fShowerFluxConstants;
88  double fSampleTime=0.;
89  double fToffset=0.;
90  std::vector<double> fBuffBox;
92  sqlite3* fdb[5];
93  double fRandomXZShift=0.;
94  };
95 }
96 
97 namespace evgen{
98 
100  : fProjectToHeight(p.get< double >("ProjectToHeight",0.)),
101  fShowerInputFiles(p.get< std::vector< std::string > >("ShowerInputFiles")),
102  fShowerFluxConstants(p.get< std::vector< double > >("ShowerFluxConstants")),
103  fSampleTime(p.get< double >("SampleTime",0.)),
104  fToffset(p.get< double >("TimeOffset",0.)),
105  fBuffBox(p.get< std::vector< double > >("BufferBox",{0.0, 0.0, 0.0, 0.0, 0.0, 0.0})),
106  fShowerAreaExtension(p.get< double >("ShowerAreaExtension",0.)),
107  fRandomXZShift(p.get< double >("RandomXZShift",0.))
108  {
109 
110  if(fShowerInputFiles.size() != fShowerFluxConstants.size() || fShowerInputFiles.size()==0 || fShowerFluxConstants.size()==0)
111  throw cet::exception("CORSIKAGen") << "ShowerInputFiles and ShowerFluxConstants have different or invalid sizes!"<<"\n";
112  fShowerInputs=fShowerInputFiles.size();
113 
114  if(fSampleTime==0.) throw cet::exception("CORSIKAGen") << "SampleTime not set!";
115 
116  if(fProjectToHeight==0.) mf::LogInfo("CORSIKAGen")<<"Using 0. for fProjectToHeight!"
117  ;
118  // create a default random engine; obtain the random seed from NuRandomService,
119  // unless overridden in configuration with key "Seed"
120  art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this, "HepJamesRandom", "gen", p, { "Seed", "SeedGenerator" });
121  art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this, "HepJamesRandom", "pois", p, "SeedPoisson");
122 
123  this->reconfigure(p);
124 
125  this->openDBs();
126  this->populateNShowers();
127  this->populateTOffset();
128 
129  produces< std::vector<simb::MCTruth> >();
130  produces< sumdata::RunData, art::InRun >();
131 
132  }
133 
135  for(int i=0; i<fShowerInputs; i++){
136  sqlite3_close(fdb[i]);
137  }
138  //cleanup temp files
139  fIFDH->cleanup();
140  }
141 
142  void CORSIKAGen::ProjectToBoxEdge( const double xyz[],
143  const double indxyz[],
144  const double xlo,
145  const double xhi,
146  const double ylo,
147  const double yhi,
148  const double zlo,
149  const double zhi,
150  double xyzout[] ){
151 
152 
153  //we want to project backwards, so take mirror of momentum
154  const double dxyz[3]={-indxyz[0],-indxyz[1],-indxyz[2]};
155 
156  // Compute the distances to the x/y/z walls
157  double dx = 99.E99;
158  double dy = 99.E99;
159  double dz = 99.E99;
160  if (dxyz[0] > 0.0) { dx = (xhi-xyz[0])/dxyz[0]; }
161  else if (dxyz[0] < 0.0) { dx = (xlo-xyz[0])/dxyz[0]; }
162  if (dxyz[1] > 0.0) { dy = (yhi-xyz[1])/dxyz[1]; }
163  else if (dxyz[1] < 0.0) { dy = (ylo-xyz[1])/dxyz[1]; }
164  if (dxyz[2] > 0.0) { dz = (zhi-xyz[2])/dxyz[2]; }
165  else if (dxyz[2] < 0.0) { dz = (zlo-xyz[2])/dxyz[2]; }
166 
167 
168  // Choose the shortest distance
169  double d = 0.0;
170  if (dx < dy && dx < dz) d = dx;
171  else if (dy < dz && dy < dx) d = dy;
172  else if (dz < dx && dz < dy) d = dz;
173 
174  // Make the step
175  for (int i = 0; i < 3; ++i) {
176  xyzout[i] = xyz[i] + dxyz[i]*d;
177  }
178 
179  }
180 
182  //choose files based on fShowerInputFiles, copy them with ifdh, open them
183  // for c2: statement is unused
184  //sqlite3_stmt *statement;
185  //get rng engine
187  CLHEP::HepRandomEngine &engine = rng->getEngine("gen");
188  CLHEP::RandFlat flat(engine);
189 
190  //setup ifdh object
191  if ( ! fIFDH ) fIFDH = new ifdh_ns::ifdh;
192  const char* ifdh_debug_env = std::getenv("IFDH_DEBUG_LEVEL");
193  if ( ifdh_debug_env ) {
194  mf::LogInfo("CORSIKAGen") << "IFDH_DEBUG_LEVEL: " << ifdh_debug_env<<"\n";
195  fIFDH->set_debug(ifdh_debug_env);
196  }
197 
198  //get ifdh path for each file in fShowerInputFiles, put into selectedflist
199  //if 1 file returned, use that file
200  //if >1 file returned, randomly select one file
201  //if 0 returned, make exeption for missing files
202  std::vector<std::pair<std::string,long>> selectedflist;
203  for(int i=0; i<fShowerInputs; i++){
204  if(fShowerInputFiles[i].find("*")==std::string::npos){
205  //if there are no wildcards, don't call findMatchingFiles
206  selectedflist.push_back(std::make_pair(fShowerInputFiles[i],0));
207  mf::LogInfo("CorsikaGen") << "Selected"<<selectedflist.back().first<<"\n";
208  }else{
209  //use findMatchingFiles
210  std::vector<std::pair<std::string,long>> flist;
211  std::string path(gSystem->DirName(fShowerInputFiles[i].c_str()));
212  std::string pattern(gSystem->BaseName(fShowerInputFiles[i].c_str()));
213  flist = fIFDH->findMatchingFiles(path,pattern);
214  unsigned int selIndex=-1;
215  if(flist.size()==1){ //0th element is the search path:pattern
216  selIndex=0;
217  }else if(flist.size()>1){
218  selIndex= (unsigned int) (flat()*(flist.size()-1)+0.5); //rnd with rounding, dont allow picking the 0th element
219  }else{
220  throw cet::exception("CORSIKAGen") << "No files returned for path:pattern: "<<path<<":"<<pattern<<std::endl;
221  }
222  selectedflist.push_back(flist[selIndex]);
223  mf::LogInfo("CorsikaGen") << "For "<<fShowerInputFiles[i]<<":"<<pattern
224  <<"\nFound "<< flist.size() << " candidate files"
225  <<"\nChoosing file number "<< selIndex << "\n"
226  <<"\nSelected "<<selectedflist.back().first<<"\n";
227  }
228 
229  }
230 
231  //do the fetching, store local filepaths in locallist
232  std::vector<std::string> locallist;
233  for(unsigned int i=0; i<selectedflist.size(); i++){
234  mf::LogInfo("CorsikaGen")
235  << "Fetching: "<<selectedflist[i].first<<" "<<selectedflist[i].second<<"\n";
236  std::string fetchedfile(fIFDH->fetchInput(selectedflist[i].first));
237  LOG_DEBUG("CorsikaGen") << "Fetched; local path: "<<fetchedfile;
238  locallist.push_back(fetchedfile);
239  }
240 
241  //open the files in fShowerInputFilesLocalPaths with sqlite3
242  for(unsigned int i=0; i<locallist.size(); i++){
243  //prepare and execute statement to attach db file
244  int res=sqlite3_open(locallist[i].c_str(),&fdb[i]);
245  if (res!= SQLITE_OK)
246  throw cet::exception("CORSIKAGen") << "Error opening db: (" <<locallist[i].c_str()<<") ("<<res<<"): " << sqlite3_errmsg(fdb[i]) << "; memory used:<<"<<sqlite3_memory_used()<<"/"<<sqlite3_memory_highwater(0)<<"\n";
247  else
248  mf::LogInfo("CORSIKAGen")<<"Attached db "<< locallist[i]<<"\n";
249  }
250  }
251 
252  double CORSIKAGen::wrapvar( const double var, const double low, const double high){
253  //wrap variable so that it's always between low and high
254  return (var - (high - low) * floor(var/(high-low))) + low;
255  }
256 
257  double CORSIKAGen::wrapvarBoxNo( const double var, const double low, const double high, int& boxno){
258  //wrap variable so that it's always between low and high
259  boxno=int(floor(var/(high-low)));
260  return (var - (high - low) * floor(var/(high-low))) + low;
261  }
262 
264  //populate TOffset_corsika by finding minimum ParticleTime from db file
265 
266  sqlite3_stmt *statement;
267  const std::string kStatement("select min(t) from particles");
268  double t=0.;
269 
270  for(int i=0; i<fShowerInputs; i++){
271  //build and do query to get run min(t) from each db
272  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
273  int res=0;
274  res = sqlite3_step(statement);
275  if ( res == SQLITE_ROW ){
276  t=sqlite3_column_double(statement,0);
277  mf::LogInfo("CORSIKAGen")<<"For showers input "<< i<<" found particles.min(t)="<<t<<"\n";
278  if (i==0 || t<fToffset_corsika) fToffset_corsika=t;
279  }else{
280  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
281  }
282  }else{
283  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
284  }
285  }
286 
287  mf::LogInfo("CORSIKAGen")<<"Found corsika timeoffset [ns]: "<< fToffset_corsika<<"\n";
288  }
289 
291  //populate vector of the number of showers per event based on:
292  //AREA the showers are being distributed over
293  //TIME of the event (fSampleTime)
294  //flux constants that determine the overall normalizations (fShowerFluxConstants)
295  //Energy range over which the sample was generated (ERANGE_*)
296  //power spectrum over which the sample was generated (ESLOPE)
297 
298 
299  //compute shower area based on the maximal x,z dimensions of cryostat boundaries + fShowerAreaExtension
301  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
302  double bounds[6] = {0.};
303  geom->CryostatBoundaries(bounds, c);
304  for (unsigned int bnd = 0; bnd<6; bnd++){
305  mf::LogVerbatim("CORSIKAGen")<<"Cryo Boundary: "<<bnd<<"="<<bounds[bnd]<<" ( + Buffer="<<fBuffBox[bnd]<<")\n";
306  if(fabs(bounds[bnd])>fabs(fShowerBounds[bnd])){
307  fShowerBounds[bnd]=bounds[bnd];
308  }
309  }
310  }
311  //add on fShowerAreaExtension without being clever
316 
317  double showersArea=(fShowerBounds[1]/100-fShowerBounds[0]/100)*(fShowerBounds[5]/100-fShowerBounds[4]/100);
318 
319  mf::LogInfo("CORSIKAGen")
320  << "Area extended by : "<<fShowerAreaExtension
321  <<"\nShowers to be distributed betweeen: x="<<fShowerBounds[0]<<","<<fShowerBounds[1]
322  <<" & z="<<fShowerBounds[4]<<","<<fShowerBounds[5]
323  <<"\nShowers to be distributed with random XZ shift: "<<fRandomXZShift<<" cm"
324  <<"\nShowers to be distributed over area: "<<showersArea<<" m^2"
325  <<"\nShowers to be distributed over time: "<<fSampleTime<<" s"
326  <<"\nShowers to be distributed with time offset: "<<fToffset<<" s"
327  <<"\nShowers to be distributed at y: "<<fShowerBounds[3]<<" cm"
328  ;
329 
330  //db variables
331  sqlite3_stmt *statement;
332  const std::string kStatement("select erange_high,erange_low,eslope,nshow from input");
333  double upperLimitOfEnergyRange=0.,lowerLimitOfEnergyRange=0.,energySlope=0.,oneMinusGamma=0.,EiToOneMinusGamma=0.,EfToOneMinusGamma=0.;
334 
335  /*
336  // get random number generator engine... currently unused here
337  art::ServiceHandle<art::RandomNumberGenerator> rng;
338  CLHEP::HepRandomEngine &engine = rng->getEngine("gen");
339  CLHEP::RandFlat flat(engine);
340  */
341 
342  for(int i=0; i<fShowerInputs; i++){
343  //build and do query to get run info from databases
344  // double thisrnd=flat();//need a new random number for each query
345  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
346  int res=0;
347  res = sqlite3_step(statement);
348  if ( res == SQLITE_ROW ){
349  upperLimitOfEnergyRange=sqlite3_column_double(statement,0);
350  lowerLimitOfEnergyRange=sqlite3_column_double(statement,1);
351  energySlope = sqlite3_column_double(statement,2);
352  fMaxShowers.push_back(sqlite3_column_int(statement,3));
353  oneMinusGamma = 1 + energySlope;
354  EiToOneMinusGamma = pow(lowerLimitOfEnergyRange, oneMinusGamma);
355  EfToOneMinusGamma = pow(upperLimitOfEnergyRange, oneMinusGamma);
356  mf::LogVerbatim("CORSIKAGen")<<"For showers input "<< i<<" found e_hi="<<upperLimitOfEnergyRange<<", e_lo="<<lowerLimitOfEnergyRange<<", slope="<<energySlope<<", k="<<fShowerFluxConstants[i]<<"\n";
357  }else{
358  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
359  }
360  }else{
361  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
362  }
363 
364  //this is computed, how?
365  double NShowers=( M_PI * showersArea * fShowerFluxConstants[i] * (EfToOneMinusGamma - EiToOneMinusGamma) / oneMinusGamma )*fSampleTime;
366  fNShowersPerEvent.push_back(NShowers);
367  mf::LogVerbatim("CORSIKAGen")<<"For showers input "<< i
368  <<" the number of showers per event is "<<(int)NShowers<<"\n";
369  }
370  }
371 
373  //for each input, randomly pull fNShowersPerEvent[i] showers from the Particles table
374  //and randomly place them in time (between -fSampleTime/2 and fSampleTime/2)
375  //wrap their positions based on the size of the area under consideration
376  //based on http://nusoft.fnal.gov/larsoft/doxsvn/html/CRYHelper_8cxx_source.html (Sample)
377 
378  //query from sqlite db with select * from particles where shower in (select id from showers ORDER BY substr(id*0.51123124141,length(id)+2) limit 100000) ORDER BY substr(shower*0.51123124141,length(shower)+2);
379  //where 0.51123124141 is a random seed to allow randomly selecting rows and should be randomly generated for each query
380  //the inner order by is to select randomly from the possible shower id's
381  //the outer order by is to make sure the shower numbers are ordered randomly (without this, the showers always come out ordered by shower number
382  //and 100000 is the number of showers to be selected at random and needs to be less than the number of showers in the showers table
383 
384  //TDatabasePDG is for looking up particle masses
385  static TDatabasePDG* pdgt = TDatabasePDG::Instance();
386 
387  //db variables
388  sqlite3_stmt *statement;
389  const TString kStatement("select shower,pdg,px,py,pz,x,z,t,e from particles where shower in (select id from showers ORDER BY substr(id*%f,length(id)+2) limit %d) ORDER BY substr(shower*%f,length(shower)+2)");
390 
391  //get rng engine
393  CLHEP::HepRandomEngine &engine = rng->getEngine("gen");
394  CLHEP::RandFlat flat(engine);
395 
396  CLHEP::HepRandomEngine &engine_pois = rng->getEngine("pois");
397  CLHEP::RandPoissonQ randpois(engine_pois);
398 
399  // get geometry and figure where to project particles to, based on CRYHelper
401  double x1, x2;
402  double y1, y2;
403  double z1, z2;
404  geom->WorldBox(&x1, &x2, &y1, &y2, &z1, &z2);
405 
406  // make the world box slightly smaller so that the projection to
407  // the edge avoids possible rounding errors later on with Geant4
408  double fBoxDelta=1.e-5;
409  x1 += fBoxDelta;
410  x2 -= fBoxDelta;
411  y1 += fBoxDelta;
412  y2 = fProjectToHeight;
413  z1 += fBoxDelta;
414  z2 -= fBoxDelta;
415 
416  //populate mctruth
417  int ntotalCtr=0; //count number of particles added to mctruth
418  int lastShower=0; //keep track of last shower id so that t can be randomized on every new shower
419  int nShowerCntr=0; //keep track of how many showers are left to be added to mctruth
420  int nShowerQry=0; //number of showers to query from db
421  int shower,pdg;
422  double px,py,pz,x,z,tParticleTime,etot,showerTime=0.,showerTimex=0.,showerTimez=0.,showerXOffset=0.,showerZOffset=0.,t;
423  for(int i=0; i<fShowerInputs; i++){
424  nShowerCntr=randpois.fire(fNShowersPerEvent[i]);
425  mf::LogInfo("CORSIKAGEN") << " Shower input " << i << " with mean " << fNShowersPerEvent[i] << " generating " << nShowerCntr;
426 
427  while(nShowerCntr>0){
428  //how many showers should we query?
429  if(nShowerCntr>fMaxShowers[i]){
430  nShowerQry=fMaxShowers[i]; //take the group size
431  }else{
432  nShowerQry=nShowerCntr; //take the rest that are needed
433  }
434  //build and do query to get nshowers
435  double thisrnd=flat(); //need a new random number for each query
436  TString kthisStatement=TString::Format(kStatement.Data(),thisrnd,nShowerQry,thisrnd);
437  LOG_DEBUG("CORSIKAGen")<<"Executing: "<<kthisStatement;
438  if ( sqlite3_prepare(fdb[i], kthisStatement.Data(), -1, &statement, 0 ) == SQLITE_OK ){
439  int res=0;
440  //loop over database rows, pushing particles into mctruth object
441  while(1){
442  res = sqlite3_step(statement);
443  if ( res == SQLITE_ROW ){
444  shower=sqlite3_column_int(statement,0);
445  if(shower!=lastShower){
446  //each new shower gets its own random time and position offsets
447  showerTime=1e9*(flat()*fSampleTime); //converting from s to ns
448  showerTimex=1e9*(flat()*fSampleTime); //converting from s to ns
449  showerTimez=1e9*(flat()*fSampleTime); //converting from s to ns
450  //and a random offset in both z and x controlled by the fRandomXZShift parameter
451  showerXOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
452  showerZOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
453  }
454  pdg=sqlite3_column_int(statement,1);
455  //get mass for this particle
456  double m = 0.; // in GeV
457  TParticlePDG* pdgp = pdgt->GetParticle(pdg);
458  if (pdgp) m = pdgp->Mass();
459 
460  //Note: position/momentum in db have north=-x and west=+z, rotate so that +z is north and +x is west
461  //get momentum components
462  px=sqlite3_column_double(statement,4);//uboone x=Particlez
463  py=sqlite3_column_double(statement,3);
464  pz=-sqlite3_column_double(statement,2);//uboone z=-Particlex
465  etot=sqlite3_column_double(statement,8);
466 
467  //get/calculate position components
468  int boxnoX=0,boxnoZ=0;
469  x=wrapvarBoxNo(sqlite3_column_double(statement,6)+showerXOffset,fShowerBounds[0],fShowerBounds[1],boxnoX);
470  z=wrapvarBoxNo(-sqlite3_column_double(statement,5)+showerZOffset,fShowerBounds[4],fShowerBounds[5],boxnoZ);
471  tParticleTime=sqlite3_column_double(statement,7); //time offset, includes propagation time from top of atmosphere
472  //actual particle time is particle surface arrival time
473  //+ shower start time
474  //+ global offset (fcl parameter, in s)
475  //- propagation time through atmosphere
476  //+ boxNo{X,Z} time offset to make grid boxes have different shower times
477  t=tParticleTime+showerTime+(1e9*fToffset)-fToffset_corsika + showerTimex*boxnoX + showerTimez*boxnoZ;
478  //wrap surface arrival so that it's in the desired time window
479  t=wrapvar(t,(1e9*fToffset),1e9*(fToffset+fSampleTime));
480 
481  simb::MCParticle p(ntotalCtr,pdg,"primary",-200,m,1);
482 
483  //project back to wordvol/fProjectToHeight
484  double xyzo[3];
485  double x0[3]={x,fShowerBounds[3],z};
486  double dx[3]={px,py,pz};
487  this->ProjectToBoxEdge(x0, dx, x1, x2, y1, y2, z1, z2, xyzo);
488 
489  TLorentzVector pos(xyzo[0],xyzo[1],xyzo[2],t);// time needs to be in ns to match GENIE, etc
490  TLorentzVector mom(px,py,pz,etot);
491  p.AddTrajectoryPoint(pos,mom);
492  mctruth.Add(p);
493  ntotalCtr++;
494  lastShower=shower;
495  }else if ( res == SQLITE_DONE ){
496  break;
497  }else{
498  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
499  }
500  }
501  }else{
502  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kthisStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
503  }
504  nShowerCntr=nShowerCntr-nShowerQry;
505  }
506  }
507  }
508 
510 
511  return;
512  }
513 
516 
517 
518  }
519 
521 
522  // grab the geometry object to see what geometry we are using
524 
525  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
526 
527  run.put(std::move(runcol));
528 
529  return;
530  }
531 
533  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
534 
536 
537  simb::MCTruth truth;
539 
540  simb::MCTruth pretruth;
541  GetSample(pretruth);
542  mf::LogInfo("CORSIKAGen")<<"GetSample number of particles returned: "<<pretruth.NParticles()<<"\n";
543  // loop over particles in the truth object
544  for(int i = 0; i < pretruth.NParticles(); ++i){
545  simb::MCParticle particle = pretruth.GetParticle(i);
546  const TLorentzVector& v4 = particle.Position();
547  const TLorentzVector& p4 = particle.Momentum();
548  double x0[3] = {v4.X(), v4.Y(), v4.Z() };
549  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
550 
551  // now check if the particle goes through any cryostat in the detector
552  // if so, add it to the truth object.
553  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
554  double bounds[6] = {0.};
555  geom->CryostatBoundaries(bounds, c);
556 
557  //add a buffer box around the cryostat bounds to increase the acceptance and account for scattering
558  //By default, the buffer box has zero size
559  for (unsigned int cb=0; cb<6; cb++)
560  bounds[cb] = bounds[cb]+fBuffBox[cb];
561 
562  //calculate the intersection point with each cryostat surface
563  bool intersects_cryo = false;
564  for (int bnd=0; bnd!=6; ++bnd) {
565  if (bnd<2) {
566  double p2[3] = {bounds[bnd], x0[1] + (dx[1]/dx[0])*(bounds[bnd] - x0[0]), x0[2] + (dx[2]/dx[0])*(bounds[bnd] - x0[0])};
567  if ( p2[1] >= bounds[2] && p2[1] <= bounds[3] &&
568  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
569  intersects_cryo = true;
570  break;
571  }
572  }
573  else if (bnd>=2 && bnd<4) {
574  double p2[3] = {x0[0] + (dx[0]/dx[1])*(bounds[bnd] - x0[1]), bounds[bnd], x0[2] + (dx[2]/dx[1])*(bounds[bnd] - x0[1])};
575  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
576  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
577  intersects_cryo = true;
578  break;
579  }
580  }
581  else if (bnd>=4) {
582  double p2[3] = {x0[0] + (dx[0]/dx[2])*(bounds[bnd] - x0[2]), x0[1] + (dx[1]/dx[2])*(bounds[bnd] - x0[2]), bounds[bnd]};
583  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
584  p2[1] >= bounds[2] && p2[1] <= bounds[3] ) {
585  intersects_cryo = true;
586  break;
587  }
588  }
589  }
590 
591  if (intersects_cryo){
592  truth.Add(particle);
593  break; //leave loop over cryostats to avoid adding particle multiple times
594  }// end if particle goes into a cryostat
595  }// end loop over cryostats in the detector
596 
597  }// loop on particles
598 
599  mf::LogInfo("CORSIKAGen")<<"Number of particles from getsample crossing cryostat + bounding box: "<<truth.NParticles()<<"\n";
600 
601  truthcol->push_back(truth);
602  evt.put(std::move(truthcol));
603 
604  return;
605  }// end produce
606 
607 }// end namespace
608 
609 
610 namespace evgen{
611 
613 
614 }
615 
616 #endif // EVGEN_CORSIKAGen_H
617 
Float_t x
Definition: compare.C:6
double fToffset
Time offset of sample, defaults to zero (no offset) [s].
double wrapvar(const double var, const double low, const double high)
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:223
double fSampleTime
Duration of sample [s].
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:261
double fProjectToHeight
Height to which particles will be projected [cm].
Encapsulate the construction of a single cyostat.
Float_t y1[n_points_granero]
Definition: compare.C:5
Collect all the geometry header files together.
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:78
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
void WorldBox(double *xlo, double *xhi, double *ylo, double *yhi, double *zlo, double *zhi) const
Fills the arguments with the boundaries of the world.
Float_t x1[n_points_granero]
Definition: compare.C:5
Double_t z
Definition: plot.C:279
A module to check the results from the Monte Carlo generator.
STL namespace.
int NParticles() const
Definition: MCTruth.h:72
ifdh_ns::ifdh * fIFDH
(optional) flux file handling
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Particle class.
void Add(simb::MCParticle &part)
Definition: MCTruth.h:77
std::vector< std::string > fShowerInputFiles
Set of CORSIKA shower data files to use.
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
Definition: Run.h:30
Float_t y2[n_points_geant4]
Definition: compare.C:26
int fShowerInputs
Number of shower inputs to process from.
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
base_engine_t & getEngine() const
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
std::vector< double > fShowerFluxConstants
Set of flux constants to be associated with each shower data file.
void ProjectToBoxEdge(const double xyz[], const double dxyz[], const double xlo, const double xhi, const double ylo, const double yhi, const double zlo, const double zhi, double xyzout[])
Float_t d
Definition: plot.C:237
double wrapvarBoxNo(const double var, const double low, const double high, int &boxno)
std::vector< int > fMaxShowers
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:73
double fToffset_corsika
Timing offset to account for propagation time through atmosphere, populated from db.
void produce(art::Event &evt)
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:224
double fRandomXZShift
Each shower will be shifted by a random amount in xz so that showers won&#39;t repeatedly sample the same...
#define LOG_DEBUG(id)
double fShowerBounds[6]
Boundaries of area over which showers are to be distributed.
CORSIKAGen(fhicl::ParameterSet const &pset)
void beginRun(art::Run &run)
TCEvent evt
Definition: DataStructs.cxx:5
Event generator information.
Definition: MCTruth.h:30
Float_t x2[n_points_geant4]
Definition: compare.C:26
void GetSample(simb::MCTruth &)
std::vector< double > fNShowersPerEvent
Number of showers to put in each event of duration fSampleTime; one per showerinput.
Namespace collecting geometry-related classes utilities.
void reconfigure(fhicl::ParameterSet const &p)
Event Generation using GENIE, cosmics or single particles.
BoundingBox bounds(int x, int y, int w, int h)
Definition: main.cpp:37
double fShowerAreaExtension
Extend distribution of corsika particles in x,z by this much (e.g. 1000 will extend 10 m in -x...
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Cosmic rays.
Definition: MCTruth.h:22
std::vector< double > fBuffBox
Buffer box extensions to cryostat in each direction (6 of them: x_lo,x_hi,y_lo,y_hi,z_lo,z_hi) [cm].