CORSIKAGendp_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file CORSIKAGendp_module.cc
3 /// \brief Generator for cosmic-ray secondaries based on pre-generated CORSIKA shower databases.
4 /// \Basic a copy of larsim/EventGenerator/CORSIKA/CORSIKAGen_module.cc by mike bass adapted for Dual-Phase purposes
5 /// \author ascarpel@apc.in2p3.fr
6 ////////////////////////////////////////////////////////////////////////
7 //#ifndef EVGEN_CORSIKAGendp_H
8 //#define EVGEN_CORSIKAGendp_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"
24 #include "art_root_io/TFileService.h"
25 #include "art_root_io/TFileDirectory.h"
27 
28 // art extensions
29 #include "nurandom/RandomUtils/NuRandomService.h"
30 
31 // larsoft includes
34 #include "nugen/EventGeneratorBase/evgenbase.h"
39 
40 #include <sqlite3.h>
41 #include "CLHEP/Random/RandFlat.h"
42 #include "CLHEP/Random/RandPoissonQ.h"
43 #include "ifdh.h" //to handle flux files
44 
45 // c/c++ includes
46 #include <sqlite3.h>
47 #include <memory>
48 #include <stdio.h>
49 #include <dirent.h>
50 
51 namespace evgendp{
52 
53  /// A module to check the results from the Monte Carlo generator
54  class CORSIKAGendp : public art::EDProducer {
55  public:
56  explicit CORSIKAGendp(fhicl::ParameterSet const& pset);
57  virtual ~CORSIKAGendp();
58 
59  void produce(art::Event& evt);
60  void beginJob();
61  void beginRun(art::Run& run);
62  void reconfigure(fhicl::ParameterSet const& p);
63 
64 
65  private:
66  void openDBs(std::string const& module_label);
67  void populateNShowers();
68  void populateTOffset();
69  void GetSample(simb::MCTruth&);
70  double wrapvar( const double var, const double low, const double high);
71  double wrapvarBoxNo( const double var, const double low, const double high, int& boxno);
72  void DoRotation(double vec[]);
73  void ProjectToBoxEdge(const double xyz[],
74  const double dxyz[],
75  const double xlo,
76  const double xhi,
77  const double ylo,
78  const double yhi,
79  const double zlo,
80  const double zhi,
81  double xyzout[]);
82 
83  int fShowerInputs=0; ///< Number of shower inputs to process from
84  std::vector<double> fNShowersPerEvent; ///< Number of showers to put in each event of duration fSampleTime; one per showerinput
85  std::vector<int> fMaxShowers; //< Max number of showers to query, one per showerinput
86  double fShowerBounds[6]={0.,0.,0.,0.,0.,0.}; ///< Boundaries of area over which showers are to be distributed
87  double fToffset_corsika=0.; ///< Timing offset to account for propagation time through atmosphere, populated from db
88  ifdh_ns::ifdh* fIFDH=0; ///< (optional) flux file handling
89 
90  //fcl parameters
91  double fProjectToHeight=0.; ///< Height to which particles will be projected [cm]
92  std::vector< std::string > fShowerInputFiles; ///< Set of CORSIKA shower data files to use
93  std::vector< double > fShowerFluxConstants; ///< Set of flux constants to be associated with each shower data file
94  double fSampleTime=0.; ///< Duration of sample [s]
95  double fToffset=0.; ///< Time offset of sample, defaults to zero (no offset) [s]
96  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]
97  double fShowerAreaExtension=0.; ///< Extend distribution of corsika particles in x,z by this much (e.g. 1000 will extend 10 m in -x, +x, -z, and +z) [cm]
98  sqlite3* fdb[5]; ///< Pointers to sqlite3 database object, max of 5
99  double fRandomXZShift=0.; ///< Each shower will be shifted by a random amount in xz so that showers won't repeatedly sample the same space [cm]
100  bool fDoRotation; //do rotation?
101  bool fUseIFDH; ///<< use ifdh protocol?
102 
103  CLHEP::HepRandomEngine& fGenEngine;
104  CLHEP::HepRandomEngine& fPoisEngine;
105  };
106 }
107 
108 namespace evgendp{
109 
111  : EDProducer{p},
112  fProjectToHeight(p.get< double >("ProjectToHeight",0.)),
113  fShowerInputFiles(p.get< std::vector< std::string > >("ShowerInputFiles")),
114  fShowerFluxConstants(p.get< std::vector< double > >("ShowerFluxConstants")),
115  fSampleTime(p.get< double >("SampleTime",0.)),
116  fToffset(p.get< double >("TimeOffset",0.)),
117  fBuffBox(p.get< std::vector< double > >("BufferBox",{0.0, 0.0, 0.0, 0.0, 0.0, 0.0})),
118  fShowerAreaExtension(p.get< double >("ShowerAreaExtension",0.)),
119  fRandomXZShift(p.get< double >("RandomXZShift",0.)),
120  fDoRotation(p.get< bool >("DoRotation")),
121  fUseIFDH(p.get< bool >("UseIFDH")),
122  // create a default random engine; obtain the random seed from NuRandomService,
123  // unless overridden in configuration with key "Seed"
125  ->createEngine(*this, "HepJamesRandom", "gen", p, { "Seed", "SeedGenerator" })),
127  ->createEngine(*this, "HepJamesRandom", "pois", p, "SeedPoisson"))
128  {
129 
130  if(fShowerInputFiles.size() != fShowerFluxConstants.size() || fShowerInputFiles.size()==0 || fShowerFluxConstants.size()==0)
131  throw cet::exception("CORSIKAGendp") << "ShowerInputFiles and ShowerFluxConstants have different or invalid sizes!"<<"\n";
133 
134  if(fSampleTime==0.) throw cet::exception("CORSIKAGendp") << "SampleTime not set!";
135 
136  if(fProjectToHeight==0.) mf::LogInfo("CORSIKAGendp")<<"Using 0. for fProjectToHeight!"
137  ;
138 
139  this->reconfigure(p);
140 
141  this->openDBs(p.get<std::string>("module_label"));
142  this->populateNShowers();
143  this->populateTOffset();
144 
145  produces< std::vector<simb::MCTruth> >();
146  produces< sumdata::RunData, art::InRun >();
147 
148  }
149 
151  for(int i=0; i<fShowerInputs; i++){
152  sqlite3_close(fdb[i]);
153  }
154  //cleanup temp files
155  fIFDH->cleanup();
156  }
157 
158  void CORSIKAGendp::DoRotation(double vec[] ){
159  //rotate a vector of 90^ around z (ProtoDUNE s.o.r.)
160  double dummyvec[3];
161  dummyvec[0] = vec[0];
162  dummyvec[1] = vec[1];
163  dummyvec[2] = vec[3];
164  vec[0] = dummyvec[1]; // x' = y
165  vec[1] = -dummyvec[0]; //y' = -x
166  vec[2] = dummyvec[2]; //z' = z
167 
168  return;
169  }
170 
171  void CORSIKAGendp::ProjectToBoxEdge( const double xyz[],
172  const double indxyz[],
173  const double xlo,
174  const double xhi,
175  const double ylo,
176  const double yhi,
177  const double zlo,
178  const double zhi,
179  double xyzout[] ){
180 
181 
182  //we want to project backwards, so take mirror of momentum
183  const double dxyz[3]={-indxyz[0],-indxyz[1],-indxyz[2]};
184 
185  // Compute the distances to the x/y/z walls
186  double dx = 99.E99;
187  double dy = 99.E99;
188  double dz = 99.E99;
189  if (dxyz[0] > 0.0) { dx = (xhi-xyz[0])/dxyz[0]; }
190  else if (dxyz[0] < 0.0) { dx = (xlo-xyz[0])/dxyz[0]; }
191  if (dxyz[1] > 0.0) { dy = (yhi-xyz[1])/dxyz[1]; }
192  else if (dxyz[1] < 0.0) { dy = (ylo-xyz[1])/dxyz[1]; }
193  if (dxyz[2] > 0.0) { dz = (zhi-xyz[2])/dxyz[2]; }
194  else if (dxyz[2] < 0.0) { dz = (zlo-xyz[2])/dxyz[2]; }
195 
196 
197  // Choose the shortest distance
198  double d = 0.0;
199  if (dx < dy && dx < dz) d = dx;
200  else if (dy < dz && dy < dx) d = dy;
201  else if (dz < dx && dz < dy) d = dz;
202 
203  // Make the step
204  for (int i = 0; i < 3; ++i) {
205  xyzout[i] = xyz[i] + dxyz[i]*d;
206  }
207 
208  }
209 
210  void CORSIKAGendp::openDBs(std::string const& module_label){
211  //choose files based on fShowerInputFiles, copy them with ifdh, open them
212  //sqlite3_stmt *statement;
213  CLHEP::RandFlat flat(fGenEngine);
214 
215  //setup ifdh object
216  if(fUseIFDH){
217  if ( ! fIFDH ) fIFDH = new ifdh_ns::ifdh;
218  const char* ifdh_debug_env = std::getenv("IFDH_DEBUG_LEVEL");
219  if ( ifdh_debug_env ) {
220  mf::LogInfo("CORSIKAGendp") << "IFDH_DEBUG_LEVEL: " << ifdh_debug_env<<"\n";
221  fIFDH->set_debug(ifdh_debug_env);
222  }
223  }
224 
225  //get ifdh path for each file in fShowerInputFiles, put into selectedflist
226  //if 1 file returned, use that file
227  //if >1 file returned, randomly select one file
228  //if 0 returned, make exeption for missing files
229  std::vector<std::pair<std::string,long>> selectedflist;
230  for(int i=0; i<fShowerInputs; i++)
231  {
232  if(fShowerInputFiles[i].find("*")==std::string::npos)
233  {
234  //if there are no wildcards, don't call findMatchingFiles
235  selectedflist.push_back(std::make_pair(fShowerInputFiles[i],0));
236  mf::LogInfo("CorsikaGendp") << "Selected"<<selectedflist.back().first<<"\n";
237  }
238  else
239  {
240  //use findMatchingFiles
241  std::vector<std::pair<std::string,long>> flist;
242  std::string path(gSystem->DirName(fShowerInputFiles[i].c_str()));
243  std::string pattern(gSystem->BaseName(fShowerInputFiles[i].c_str()));
244  //pattern="/"+pattern;
245  std::cout << path << " " << pattern << std::endl;
246 
247  if(fUseIFDH)
248  {
249  //acess the files using ifdh
250  flist = fIFDH->findMatchingFiles(path,pattern);
251  }
252  else
253  {
254  //access the files using drent
255  auto wildcardPosition = pattern.find("*");
256  pattern = pattern.substr( 0, wildcardPosition );
257 
258  DIR *dir;
259  struct dirent *ent;
260  int index=0;
261  if ((dir = opendir( path.c_str() )) != NULL)
262  {
263  while ((ent = readdir (dir)) != NULL)
264  {
265  index++;
266  std::pair<std::string,long> name;
267  std::string parsename(ent->d_name);
268 
269  if( parsename.substr(0, wildcardPosition) == pattern )
270  {
271  name.first= path+"/"+parsename;
272  name.second = index;
273  flist.push_back( name );
274  }
275  }
276  closedir (dir);
277  }
278  else
279  {
280  throw cet::exception("Gen311") << "Can't open directory with pattern: "<<path<<":"<<pattern<<std::endl;
281  }
282  }
283 
284  unsigned int selIndex=-1;
285  if(flist.size()==1) //0th element is the search path:pattern
286  {
287  selIndex=0;
288  }
289  else if(flist.size()>1)
290  {
291  selIndex= (unsigned int) (flat()*(flist.size()-1)+0.5); //rnd with rounding, dont allow picking the 0th element
292  }
293  else
294  {
295  throw cet::exception("CORSIKAGendp") << "No files returned for path:pattern: "<<path<<":"<<pattern<<std::endl;
296  }
297 
298  selectedflist.push_back(flist[selIndex]);
299  mf::LogInfo("CorsikaGendp") << "For "<<fShowerInputFiles[i]<<":"<<pattern
300  <<"\nFound "<< flist.size() << " candidate files"
301  <<"\nChoosing file number "<< selIndex << "\n"
302  <<"\nSelected "<<selectedflist.back().first<<"\n";
303  }
304  }
305 
306  //open the files in fShowerInputFilesLocalPaths with sqlite3
307  std::vector<std::string> locallist;
308  for(unsigned int i=0; i<selectedflist.size(); i++){
309  mf::LogInfo("Gen311")
310  << "Fetching: "<<selectedflist[i].first<<" "<<selectedflist[i].second<<"\n";
311  std::string fetchedfile(selectedflist[i].first);
312  MF_LOG_DEBUG("Gen311") << "Fetched; local path: "<<fetchedfile;
313  locallist.push_back(fetchedfile);
314  }
315 
316  for(unsigned int i=0; i<locallist.size(); i++){
317  //prepare and execute statement to attach db file
318  int res=sqlite3_open(locallist[i].c_str(),&fdb[i]);
319  if (res!= SQLITE_OK)
320  throw cet::exception("Gen311") << "Error opening db: (" <<locallist[i].c_str()<<") ("<<res<<"): " << sqlite3_errmsg(fdb[i]) << "; memory used:<<"<<sqlite3_memory_used()<<"/"<<sqlite3_memory_highwater(0)<<"\n";
321  else
322  mf::LogInfo("Gen311")<<"Attached db "<< locallist[i]<<"\n";
323  }
324  }//end openDBs
325 
326  double CORSIKAGendp::wrapvar( const double var, const double low, const double high){
327  //wrap variable so that it's always between low and high
328  return (var - (high - low) * floor(var/(high-low))) + low;
329  }
330 
331  double CORSIKAGendp::wrapvarBoxNo( const double var, const double low, const double high, int& boxno){
332  //wrap variable so that it's always between low and high
333  boxno=int(floor(var/(high-low)));
334  return (var - (high - low) * floor(var/(high-low))) + low;
335  }
336 
338  //populate TOffset_corsika by finding minimum ParticleTime from db file
339 
340  sqlite3_stmt *statement;
341  const std::string kStatement("select min(t) from particles");
342  double t=0.;
343 
344  for(int i=0; i<fShowerInputs; i++){
345  //build and do query to get run min(t) from each db
346  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
347  int res=0;
348  res = sqlite3_step(statement);
349  if ( res == SQLITE_ROW ){
350  t=sqlite3_column_double(statement,0);
351  mf::LogInfo("CORSIKAGendp")<<"For showers input "<< i<<" found particles.min(t)="<<t<<"\n";
352  if (i==0 || t<fToffset_corsika) fToffset_corsika=t;
353  }else{
354  throw cet::exception("CORSIKAGendp") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
355  }
356  }else{
357  throw cet::exception("CORSIKAGendp") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
358  }
359  }
360 
361  mf::LogInfo("CORSIKAGendp")<<"Found corsika timeoffset [ns]: "<< fToffset_corsika<<"\n";
362  }
363 
365  //populate vector of the number of showers per event based on:
366  //AREA the showers are being distributed over
367  //TIME of the event (fSampleTime)
368  //flux constants that determine the overall normalizations (fShowerFluxConstants)
369  //Energy range over which the sample was generated (ERANGE_*)
370  //power spectrum over which the sample was generated (ESLOPE)
371 
372 
373  //compute shower area based on the maximal x,z dimensions of cryostat boundaries + fShowerAreaExtension
375  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
376  double bounds[6] = {0.};
377  geom->CryostatBoundaries(bounds, c);
378  for (unsigned int bnd = 0; bnd<6; bnd++){
379  mf::LogVerbatim("CORSIKAGendp")<<"Cryo Boundary: "<<bnd<<"="<<bounds[bnd]<<" ( + Buffer="<<fBuffBox[bnd]<<")\n";
380  if(fabs(bounds[bnd])>fabs(fShowerBounds[bnd])){
381  fShowerBounds[bnd]=bounds[bnd];
382  }
383  }
384  }
385  //add on fShowerAreaExtension without being clever
390 
391  double showersArea=(fShowerBounds[1]/100-fShowerBounds[0]/100)*(fShowerBounds[5]/100-fShowerBounds[4]/100);
392 
393  mf::LogInfo("CORSIKAGendp")
394  << "Area extended by : "<<fShowerAreaExtension
395  <<"\nShowers to be distributed betweeen: x="<<fShowerBounds[0]<<","<<fShowerBounds[1]
396  <<" & z="<<fShowerBounds[4]<<","<<fShowerBounds[5]
397  <<"\nShowers to be distributed with random XZ shift: "<<fRandomXZShift<<" cm"
398  <<"\nShowers to be distributed over area: "<<showersArea<<" m^2"
399  <<"\nShowers to be distributed over time: "<<fSampleTime<<" s"
400  <<"\nShowers to be distributed with time offset: "<<fToffset<<" s"
401  <<"\nShowers to be distributed at y: "<<fShowerBounds[3]<<" cm"
402  ;
403 
404  //db variables
405  sqlite3_stmt *statement;
406  const std::string kStatement("select erange_high,erange_low,eslope,nshow from input");
407  double upperLimitOfEnergyRange=0.,lowerLimitOfEnergyRange=0.,energySlope=0.,oneMinusGamma=0.,EiToOneMinusGamma=0.,EfToOneMinusGamma=0.;
408 
409  /*
410  CLHEP::RandFlat flat(fGenEngine);
411  */
412 
413  for(int i=0; i<fShowerInputs; i++){
414  //build and do query to get run info from databases
415  // double thisrnd=flat();//need a new random number for each query
416  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
417  int res=0;
418  res = sqlite3_step(statement);
419  if ( res == SQLITE_ROW ){
420  upperLimitOfEnergyRange=sqlite3_column_double(statement,0);
421  lowerLimitOfEnergyRange=sqlite3_column_double(statement,1);
422  energySlope = sqlite3_column_double(statement,2);
423  fMaxShowers.push_back(sqlite3_column_int(statement,3));
424  oneMinusGamma = 1 + energySlope;
425  EiToOneMinusGamma = pow(lowerLimitOfEnergyRange, oneMinusGamma);
426  EfToOneMinusGamma = pow(upperLimitOfEnergyRange, oneMinusGamma);
427  mf::LogVerbatim("CORSIKAGendp")<<"For showers input "<< i<<" found e_hi="<<upperLimitOfEnergyRange<<", e_lo="<<lowerLimitOfEnergyRange<<", slope="<<energySlope<<", k="<<fShowerFluxConstants[i]<<"\n";
428  }else{
429  throw cet::exception("CORSIKAGendp") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
430  }
431  }else{
432  throw cet::exception("CORSIKAGendp") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
433  }
434 
435  //this is computed, how?
436  double NShowers=( M_PI * showersArea * fShowerFluxConstants[i] * (EfToOneMinusGamma - EiToOneMinusGamma) / oneMinusGamma )*fSampleTime;
437  fNShowersPerEvent.push_back(NShowers);
438  mf::LogVerbatim("CORSIKAGendp")<<"For showers input "<< i
439  <<" the number of showers per event is "<<(int)NShowers<<"\n";
440  }
441  }
442 
444  //for each input, randomly pull fNShowersPerEvent[i] showers from the Particles table
445  //and randomly place them in time (between -fSampleTime/2 and fSampleTime/2)
446  //wrap their positions based on the size of the area under consideration
447  //based on http://nusoft.fnal.gov/larsoft/doxsvn/html/CRYHelper_8cxx_source.html (Sample)
448 
449  //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);
450  //where 0.51123124141 is a random seed to allow randomly selecting rows and should be randomly generated for each query
451  //the inner order by is to select randomly from the possible shower id's
452  //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
453  //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
454 
455  //TDatabasePDG is for looking up particle masses
456  static TDatabasePDG* pdgt = TDatabasePDG::Instance();
457 
458  //db variables
459  sqlite3_stmt *statement;
460  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)");
461 
462  //get rng engine
463  CLHEP::RandFlat flat(fGenEngine);
464  CLHEP::RandPoissonQ randpois(fPoisEngine);
465 
466  // get geometry and figure where to project particles to, based on CRYHelper
468  double x1, x2;
469  double y1, y2;
470  double z1, z2;
471  geom->WorldBox(&x1, &x2, &y1, &y2, &z1, &z2);
472 
473  // make the world box slightly smaller so that the projection to
474  // the edge avoids possible rounding errors later on with Geant4
475  double fBoxDelta=1.e-5;
476  x1 += fBoxDelta;
477  x2 -= fBoxDelta;
478  y1 += fBoxDelta;
479  y2 = fProjectToHeight;
480  z1 += fBoxDelta;
481  z2 -= fBoxDelta;
482 
483  if(fDoRotation)
484  {
485  geom->WorldBox(&y1, &y2, &x1, &x2, &z1, &z2);
486  y1 += fBoxDelta;
487  y2 -= fBoxDelta;
488  x1 += fBoxDelta;
489  x2 = fProjectToHeight;
490  z1 += fBoxDelta;
491  z2 -= fBoxDelta;
492  }
493 
494  //populate mctruth
495  int ntotalCtr=0; //count number of particles added to mctruth
496  int lastShower=0; //keep track of last shower id so that t can be randomized on every new shower
497  int nShowerCntr=0; //keep track of how many showers are left to be added to mctruth
498  int nShowerQry=0; //number of showers to query from db
499  int shower,pdg;
500  double px,py,pz,x=0, y=0,z,tParticleTime,etot,showerTime=0.,showerTimex=0.,showerTimez=0.,showerXOffset=0.,showerZOffset=0.,t;
501  for(int i=0; i<fShowerInputs; i++){
502  nShowerCntr=randpois.fire(fNShowersPerEvent[i]);
503  mf::LogInfo("CORSIKAgendp") << " Shower input " << i << " with mean " << fNShowersPerEvent[i] << " generating " << nShowerCntr;
504 
505  while(nShowerCntr>0){
506  //how many showers should we query?
507  if(nShowerCntr>fMaxShowers[i]){
508  nShowerQry=fMaxShowers[i]; //take the group size
509  }else{
510  nShowerQry=nShowerCntr; //take the rest that are needed
511  }
512  //build and do query to get nshowers
513  double thisrnd=flat(); //need a new random number for each query
514  TString kthisStatement=TString::Format(kStatement.Data(),thisrnd,nShowerQry,thisrnd);
515  MF_LOG_DEBUG("CORSIKAGendp")<<"Executing: "<<kthisStatement;
516  if ( sqlite3_prepare(fdb[i], kthisStatement.Data(), -1, &statement, 0 ) == SQLITE_OK ){
517  int res=0;
518  //loop over database rows, pushing particles into mctruth object
519  while(1){
520  res = sqlite3_step(statement);
521  if ( res == SQLITE_ROW ){
522  shower=sqlite3_column_int(statement,0);
523  if(shower!=lastShower){
524  //each new shower gets its own random time and position offsets
525  showerTime=1e9*(flat()*fSampleTime); //converting from s to ns
526  showerTimex=1e9*(flat()*fSampleTime); //converting from s to ns
527  showerTimez=1e9*(flat()*fSampleTime); //converting from s to ns
528  //and a random offset in both z and x controlled by the fRandomXZShift parameter
529  showerXOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
530  showerZOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
531  }
532  pdg=sqlite3_column_int(statement,1);
533  //get mass for this particle
534  double m = 0.; // in GeV
535  TParticlePDG* pdgp = pdgt->GetParticle(pdg);
536  if (pdgp) m = pdgp->Mass();
537 
538  //Note: position/momentum in db have north=-x and west=+z, rotate so that +z is north and +x is west
539  //get momentum components
540  if(fDoRotation)
541  {
542  py = sqlite3_column_double(statement,4);
543  px=sqlite3_column_double(statement,3);
544  }
545  else
546  {
547  px=sqlite3_column_double(statement,4);//uboone x=Particlez
548  py=sqlite3_column_double(statement,3);
549  }
550  pz=-sqlite3_column_double(statement,2);//uboone z=-Particlex
551  etot=sqlite3_column_double(statement,8);
552 
553  //get/calculate position components
554  int boxnoX=0,boxnoZ=0;
555  if(fDoRotation) y = wrapvarBoxNo(sqlite3_column_double(statement,6)+showerXOffset,fShowerBounds[0],fShowerBounds[1],boxnoX);
556  else x=wrapvarBoxNo(sqlite3_column_double(statement,6)+showerXOffset,fShowerBounds[0],fShowerBounds[1],boxnoX);
557  z=wrapvarBoxNo(-sqlite3_column_double(statement,5)+showerZOffset,fShowerBounds[4],fShowerBounds[5],boxnoZ);
558  tParticleTime=sqlite3_column_double(statement,7); //time offset, includes propagation time from top of atmosphere
559  //actual particle time is particle surface arrival time
560  //+ shower start time
561  //+ global offset (fcl parameter, in s)
562  //- propagation time through atmosphere
563  //+ boxNo{X,Z} time offset to make grid boxes have different shower times
564  t=tParticleTime+showerTime+(1e9*fToffset)-fToffset_corsika + showerTimex*boxnoX + showerTimez*boxnoZ;
565  //wrap surface arrival so that it's in the desired time window
566  t=wrapvar(t,(1e9*fToffset),1e9*(fToffset+fSampleTime));
567 
568  simb::MCParticle p(ntotalCtr,pdg,"primary",-200,m,1);
569 
570  //project back to wordvol/fProjectToHeight
571  double xyzo[3];
572  double x0[3];
573  if(fDoRotation){
574  x0[0]=fShowerBounds[3];
575  x0[1]=y;
576  x0[2]=z;
577  }
578  else{
579  x0[0]=x;
580  x0[1]=fShowerBounds[3];
581  x0[2]=z;
582  }
583 
584  double dx[3]={px,py,pz};
585 
586  // non pare funzionare molto bene...
587  // if(fDoRotation)
588  // {
589  //add an additional rotation to this vectors to match the DP geom: x'=y x and y'=-x (90^ rotation)
590  // this->DoRotation(x0);
591  // this->DoRotation(dx);
592  // }
593 
594  this->ProjectToBoxEdge(x0, dx, x1, x2, y1, y2, z1, z2, xyzo);
595 
596  TLorentzVector pos(xyzo[0],xyzo[1],xyzo[2],t);// time needs to be in ns to match GENIE, etc
597  TLorentzVector mom(px,py,pz,etot);
598  p.AddTrajectoryPoint(pos,mom);
599  mctruth.Add(p);
600  ntotalCtr++;
601  lastShower=shower;
602  }else if ( res == SQLITE_DONE ){
603  break;
604  }else{
605  throw cet::exception("CORSIKAGendp") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
606  }
607  }
608  }else{
609  throw cet::exception("CORSIKAGendp") << "Error preparing statement: (" <<kthisStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
610  }
611  nShowerCntr=nShowerCntr-nShowerQry;
612  }
613  }
614  }
615 
617 
618  return;
619  }
620 
623 
624 
625  }
626 
628 
629  // grab the geometry object to see what geometry we are using
631 
632  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
633 
634  run.put(std::move(runcol));
635 
636  return;
637  }
638 
640  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
641 
643 
644  simb::MCTruth truth;
646 
647  simb::MCTruth pretruth;
648  GetSample(pretruth);
649  mf::LogInfo("CORSIKAGendp")<<"GetSample number of particles returned: "<<pretruth.NParticles()<<"\n";
650  // loop over particles in the truth object
651 
652  //add a buffer box around the cryostat bounds to increase the acceptance and account for scattering
653  //By default, the buffer box has zero size
654 
655  for(int i = 0; i < pretruth.NParticles(); ++i){
656  simb::MCParticle particle = pretruth.GetParticle(i);
657  const TLorentzVector& v4 = particle.Position();
658  const TLorentzVector& p4 = particle.Momentum();
659  double x0[3] = {v4.X(), v4.Y(), v4.Z() };
660  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
661 
662  // now check if the particle goes through any cryostat in the detector
663  // if so, add it to the truth object.
664  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
665  double bounds[6] = {0.};
666  geom->CryostatBoundaries(bounds, c);
667 
668  //add a buffer box around the cryostat bounds to increase the acceptance and account for scattering
669  //By default, the buffer box has zero size
670  for (unsigned int cb=0; cb<6; cb++){
671  bounds[cb] = bounds[cb]+fBuffBox[cb];
672  }
673 
674  //calculate the intersection point with each cryostat surface
675  bool intersects_cryo = false;
676  for (int bnd=0; bnd!=6; ++bnd) {
677  if (bnd<2) {
678  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])};
679  if ( p2[1] >= bounds[2] && p2[1] <= bounds[3] &&
680  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
681  intersects_cryo = true;
682  break;
683  }
684  }
685  else if (bnd>=2 && bnd<4) {
686  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])};
687  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
688  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
689  intersects_cryo = true;
690  break;
691  }
692  }
693  else if (bnd>=4) {
694  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]};
695  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
696  p2[1] >= bounds[2] && p2[1] <= bounds[3] ) {
697  intersects_cryo = true;
698  break;
699  }
700  }
701  }
702 
703  if (intersects_cryo){
704  truth.Add(particle);
705  break; //leave loop over cryostats to avoid adding particle multiple times
706  }// end if particle goes into a cryostat
707  }// end loop over cryostats in the detector
708 
709  }// loop on particles
710 
711  mf::LogInfo("CORSIKAGendp")<<"Number of particles from getsample crossing cryostat + bounding box: "<<truth.NParticles()<<"\n";
712 
713  truthcol->push_back(truth);
714  evt.put(std::move(truthcol));
715 
716  return;
717  }// end produce
718 
719 }// end namespace
720 
721 
722 namespace evgendp{
723 
725 
726 }
727 
728 //#endif // EVGEN_CORSIKAGendp_H
729 ////////////////////////////////////////////////////////////////////////
static QCString name
Definition: declinfo.cpp:673
base_engine_t & createEngine(seed_t seed)
std::vector< double > fNShowersPerEvent
Number of showers to put in each event of duration fSampleTime; one per showerinput.
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:219
double fShowerBounds[6]
Boundaries of area over which showers are to be distributed.
Format
Definition: utils.h:7
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:257
CORSIKAGendp(fhicl::ParameterSet const &pset)
Encapsulate the construction of a single cyostat.
Collect all the geometry header files together.
std::vector< std::string > fShowerInputFiles
Set of CORSIKA shower data files to use.
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
double fShowerAreaExtension
Extend distribution of corsika particles in x,z by this much (e.g. 1000 will extend 10 m in -x...
std::string string
Definition: nybbler.cc:12
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.
double fSampleTime
Duration of sample [s].
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
constexpr T pow(T x)
Definition: pow.h:72
std::vector< int > fMaxShowers
CLHEP::HepRandomEngine & fPoisEngine
void openDBs(std::string const &module_label)
struct sqlite3_stmt sqlite3_stmt
double wrapvar(const double var, const double low, const double high)
void GetSample(simb::MCTruth &)
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[])
void DoRotation(double vec[])
int NParticles() const
Definition: MCTruth.h:75
string dir
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Particle class.
double fToffset_corsika
Timing offset to account for propagation time through atmosphere, populated from db.
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
art framework interface to geometry description
Definition: Run.h:17
int fShowerInputs
Number of shower inputs to process from.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void reconfigure(fhicl::ParameterSet const &p)
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
CLHEP::HepRandomEngine & fGenEngine
std::string getenv(std::string const &name)
Definition: getenv.cc:15
def move(depos, offset)
Definition: depos.py:107
double fToffset
Time offset of sample, defaults to zero (no offset) [s].
p
Definition: test.py:223
void produce(art::Event &evt)
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
std::vector< double > fShowerFluxConstants
Set of flux constants to be associated with each shower data file.
BoundingBox bounds(int x, int y, int w, int h)
Definition: main.cpp:37
bool fUseIFDH
< use ifdh protocol?
struct sqlite3 sqlite3
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:76
#define M_PI
Definition: includeROOT.h:54
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].
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
double wrapvarBoxNo(const double var, const double low, const double high, int &boxno)
void beginRun(art::Run &run)
int var
Definition: 018_def.c:9
double fProjectToHeight
Height to which particles will be projected [cm].
std::string pattern
Definition: regex_t.cc:35
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
#define MF_LOG_DEBUG(id)
list x
Definition: train.py:276
double fRandomXZShift
Each shower will be shifted by a random amount in xz so that showers won&#39;t repeatedly sample the same...
TCEvent evt
Definition: DataStructs.cxx:7
Event generator information.
Definition: MCTruth.h:32
A module to check the results from the Monte Carlo generator.
ifdh_ns::ifdh * fIFDH
(optional) flux file handling
LArSoft geometry interface.
Definition: ChannelGeo.h:16
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Cosmic rays.
Definition: MCTruth.h:24
QTextStream & endl(QTextStream &s)