All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
evgen::CORSIKAGen Class Reference

LArSoft interface to CORSIKA event generator. More...

Inheritance diagram for evgen::CORSIKAGen:

Public Member Functions

 CORSIKAGen (fhicl::ParameterSet const &pset)
 
virtual ~CORSIKAGen ()
 
void produce (art::Event &evt)
 
void beginRun (art::Run &run)
 

Private Member Functions

void openDBs (std::string const &module_label)
 
void populateNShowers ()
 
void populateTOffset ()
 
void GetSample (simb::MCTruth &)
 
double wrapvar (const double var, const double low, const double high)
 
double wrapvarBoxNo (const double var, const double low, const double high, int &boxno)
 
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[])
 Propagates a point back to the surface of a box. More...
 

Private Attributes

int fShowerInputs =0
 Number of shower inputs to process from. More...
 
std::vector< double > fNShowersPerEvent
 Number of showers to put in each event of duration fSampleTime; one per showerinput. More...
 
std::vector< int > fMaxShowers
 
double fShowerBounds [6] ={0.,0.,0.,0.,0.,0.}
 Boundaries of area over which showers are to be distributed (x(min), x(max), unused, y, z(min), z(max) ) More...
 
double fToffset_corsika =0.
 Timing offset to account for propagation time through atmosphere, populated from db. More...
 
ifdh_ns::ifdh * fIFDH =0
 (optional) flux file handling More...
 
double fProjectToHeight =0.
 Height to which particles will be projected [cm]. More...
 
std::vector< std::string > fShowerInputFiles
 Set of CORSIKA shower data files to use. More...
 
std::string fShowerCopyType
 
std::vector< double > fShowerFluxConstants
 Set of flux constants to be associated with each shower data file. More...
 
double fSampleTime =0.
 Duration of sample [s]. More...
 
double fToffset =0.
 Time offset of sample, defaults to zero (no offset) [s]. More...
 
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]. More...
 
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]. More...
 
sqlite3 * fdb [5]
 Pointers to sqlite3 database object, max of 5. More...
 
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]. More...
 
CLHEP::HepRandomEngine & fGenEngine
 
CLHEP::HepRandomEngine & fPoisEngine
 

Detailed Description

LArSoft interface to CORSIKA event generator.

Note
A presentation on this module by the original author is archived at: https://indico.fnal.gov/event/10893/contribution/3/material/slides

In CORSIKA jargon, a "shower" is the cascade of particles resulting from a primary cosmic ray interaction. This module creates a single simb::MCTruth object (stored as data product into a std::vector<simb::MCTruth> with a single entry) containing all the particles from cosmic ray showers crossing a surface above the detector.

The generation procedure consists of selecting showers from a database of pregenerated events, and then to adapt them to the parameters requested in the module configuration. Pregenerated showers are "observed" at a altitude set in CORSIKA configuration.

Databases need to be stored as files in SQLite3 format. Multiple file sources can be specified (ShowerInputFiles configuration parameter). From each source, one database file is selected and copied locally via IFDH. From each source, showers are extracted proportionally to the relative flux specified in the configuration (specified in ShowerFluxConstants, see normalization below). The actual number of showers per event and per source is extracted according to a Poisson distribution around the predicted average number of primary cosmic rays for that source.

Flux normalization

CORSIKA generates showers from each specific cosmic ray type $ A $ (e.g. iron, proton, etc.) according to a power law distribution $ \Phi_{A}(E) \propto E^{-\gamma_{A}} $ of the primary particle energy $ E $ [GeV]. When sampling pregenerated events, we bypass the normalization imposed by CORSIKA and gain complete control on it.

Within CORSIKAGen, for each source (usually each on a different primary cosmic ray type, e.g. iron, proton, etc.), the average number of generated showers is $ n_{A} = \pi S T k_{A} \int E^{-\gamma_{A}} dE $ with $ S $ the area of the surface the flux passes across, $ T $ the exposure time, the integral defined over the full energy range of the pregenerated showers in the source, and $ k_{A} $ a factor specified in the configuration (ShowerFluxConstants parameters). This is the flux of primary cosmic rays, not of the observed particles from their showers. Note that it depends on an area and a time interval, but it is uniform with respect to translations and constant in time.

As explained below, we consider only the secondary particles that cross an "observation" surface. After cosmic ray primary particles cross the flux surface ( $ S_{\Phi} $ above), they develop into showers of particles that spread across large areas. Limiting ourself to the observation of particles on a small surface $ S_{o} $ has two effects. We lose the part of the showers that misses that surface $ S_{o} $. Also, considering a span of time with multiple showers, we miss particles from other hypothetical showers whose primaries crossed outside the generation surface $ S_{\Phi} $ whose shower development would leak into $ S_{o} $: we did not simulate those showers at all! In terms of total flux of observed particles, under the assumption that the flux is uniform in space, choosing $ S_{\Phi} $ the same size as $ S_{o} $ makes the two effects get the same size: just as many particles from primaries from $ S_{\Phi} $ leak out of $ S_{o} $, as many particles from primaries from outside $ S_{\Phi} $ sneak in $ S_{o} $. In that case, counting all the particles from the primaries crossing a surface $ S_{\Phi} $ of area S regardless of where they land yields the right amount of cosmic ray secondary particles across an observation surface $ S_{o} $ also of area S. Practically, the particles landing outside $ S_{o} $ need to be recovered to preserve the correct normalization, which is described in the next section.

Surface coverage, position and timing

The surface we detect the particles through (let's call it $ S_{d} $) is defined by the smallest rectangle including all cryostats in the detector, and located at the height of the ceiling of the tallest cryostat. This surface can be increased by specifying a positive value for ShowerAreaExtension configuration parameter, in which case each side of the rectangle will be extended by that amount.

Showers are extracted one by one from the pregenerated samples and treated independently. Ideally, the detection surface $ S_{d} $ would be at the same exact altitude as the observation surface set in CORSIKA (called $ S_{o} $ above). In practice, we go the other way around, with the assumption that the shower observed at $ S_{d} $ would be very close to the one we actually generated at $ S_{o} $, and teleport the generated particles on $ S_{d} $. Since the cryostats may be just meters from the earth surface $ S_{o} $ lies on, this is an acceptable approximation.

All the particles of one shower are compelled within surface $ S_{d} $ as a first step. As explained when describing the "normalization", we need to keep all the shower particles, one way or the other. So, particles of the shower that fell out of $ S_{d} $ are repackaged into other showers and translated back in. This is equivalent to assume the primary cosmic rays originating such shower would happen outside our generation volume ( $ S_{\Phi} $), and their shower would then spill into $ S_{d} $. Since these repackaged showers are in principle independent and uncorrelated, they are assigned a random time different than the main shower, leveraging the assumption of constantness of the flux.

As for the azimuth, this module uses an approximation by setting north direction to match the z axis of LArSoft geometry (typically assumed to be the direction of the beam particle).

The particles so manipulated are then back-propagated from the observation surface to an absolute height defined by ProjectToHeight (although for particular combination of position and direction, the particles might be propagated back to the edge of the world, or even outside it).

As a final filter, only the particles whose straight projections cross any of the cryostats (with some buffer volume around, defined by BufferBox) are stored, while the other ones are discarded. Note that the actual interactions that particles from the generation surface undergo may deviate them enough to miss the cryostats anyway, and conversely particles that have been filtered out because shooting off the cryostats might have been subsequently deviated to actually cross them. This effect is not corrected for at this time.

The time of the showers is uniformly distributed within the configured time interval, defined by SampleTime starting from TimeOffset.

Configuration parameters

Random engines

Currently two random engines are used:

Definition at line 221 of file CORSIKAGen_module.cc.

Constructor & Destructor Documentation

evgen::CORSIKAGen::CORSIKAGen ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 291 of file CORSIKAGen_module.cc.

292  : EDProducer{p},
293  fProjectToHeight(p.get< double >("ProjectToHeight",0.)),
294  fShowerInputFiles(p.get< std::vector< std::string > >("ShowerInputFiles")),
295  fShowerCopyType(p.get< std::string >("ShowerCopyType", "IFDH")),
296  fShowerFluxConstants(p.get< std::vector< double > >("ShowerFluxConstants")),
297  fSampleTime(p.get< double >("SampleTime",0.)),
298  fToffset(p.get< double >("TimeOffset",0.)),
299  fBuffBox(p.get< std::vector< double > >("BufferBox",{0.0, 0.0, 0.0, 0.0, 0.0, 0.0})),
double fToffset
Time offset of sample, defaults to zero (no offset) [s].
double fSampleTime
Duration of sample [s].
double fProjectToHeight
Height to which particles will be projected [cm].
pdgs p
Definition: selectors.fcl:22
std::vector< std::string > fShowerInputFiles
Set of CORSIKA shower data files to use.
std::vector< double > fShowerFluxConstants
Set of flux constants to be associated with each shower data file.
std::string fShowerCopyType
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].
evgen::CORSIKAGen::~CORSIKAGen ( )
virtual

Definition at line 325 of file CORSIKAGen_module.cc.

325  {
326  for(int i=0; i<fShowerInputs; i++){
327  sqlite3_close(fdb[i]);
328  }
329  //cleanup temp files
330  fIFDH->cleanup();
331  }
ifdh_ns::ifdh * fIFDH
(optional) flux file handling
int fShowerInputs
Number of shower inputs to process from.
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.

Member Function Documentation

void evgen::CORSIKAGen::beginRun ( art::Run &  run)

Definition at line 714 of file CORSIKAGen_module.cc.

715  {
716  art::ServiceHandle<geo::Geometry const> geo;
717  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
718  }
void evgen::CORSIKAGen::GetSample ( simb::MCTruth &  mctruth)
private

Definition at line 562 of file CORSIKAGen_module.cc.

562  {
563  //for each input, randomly pull fNShowersPerEvent[i] showers from the Particles table
564  //and randomly place them in time (between -fSampleTime/2 and fSampleTime/2)
565  //wrap their positions based on the size of the area under consideration
566  //based on http://nusoft.fnal.gov/larsoft/doxsvn/html/CRYHelper_8cxx_source.html (Sample)
567 
568  //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);
569  //where 0.51123124141 is a random seed to allow randomly selecting rows and should be randomly generated for each query
570  //the inner order by is to select randomly from the possible shower id's
571  //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
572  //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
573 
574  //TDatabasePDG is for looking up particle masses
575  static TDatabasePDG* pdgt = TDatabasePDG::Instance();
576 
577  //db variables
578  sqlite3_stmt *statement;
579  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)");
580 
581  CLHEP::RandFlat flat(fGenEngine);
582  CLHEP::RandPoissonQ randpois(fPoisEngine);
583 
584  // get geometry and figure where to project particles to, based on CRYHelper
585  art::ServiceHandle<geo::Geometry const> geom;
586  double x1, x2;
587  double y1, y2;
588  double z1, z2;
589  geom->WorldBox(&x1, &x2, &y1, &y2, &z1, &z2);
590 
591  // make the world box slightly smaller so that the projection to
592  // the edge avoids possible rounding errors later on with Geant4
593  double fBoxDelta=1.e-5;
594  x1 += fBoxDelta;
595  x2 -= fBoxDelta;
596  y1 += fBoxDelta;
597  y2 = fProjectToHeight;
598  z1 += fBoxDelta;
599  z2 -= fBoxDelta;
600 
601  //populate mctruth
602  int ntotalCtr=0; //count number of particles added to mctruth
603  int lastShower=0; //keep track of last shower id so that t can be randomized on every new shower
604  int nShowerCntr=0; //keep track of how many showers are left to be added to mctruth
605  int nShowerQry=0; //number of showers to query from db
606  int shower,pdg;
607  double px,py,pz,x,z,tParticleTime,etot,showerTime=0.,showerTimex=0.,showerTimez=0.,showerXOffset=0.,showerZOffset=0.,t;
608  for(int i=0; i<fShowerInputs; i++){
609  nShowerCntr=randpois.fire(fNShowersPerEvent[i]);
610  mf::LogInfo("CORSIKAGEN") << " Shower input " << i << " with mean " << fNShowersPerEvent[i] << " generating " << nShowerCntr;
611 
612  while(nShowerCntr>0){
613  //how many showers should we query?
614  if(nShowerCntr>fMaxShowers[i]){
615  nShowerQry=fMaxShowers[i]; //take the group size
616  }else{
617  nShowerQry=nShowerCntr; //take the rest that are needed
618  }
619  //build and do query to get nshowers
620  double thisrnd=flat(); //need a new random number for each query
621  TString kthisStatement=TString::Format(kStatement.Data(),thisrnd,nShowerQry,thisrnd);
622  MF_LOG_DEBUG("CORSIKAGen")<<"Executing: "<<kthisStatement;
623  if ( sqlite3_prepare(fdb[i], kthisStatement.Data(), -1, &statement, 0 ) == SQLITE_OK ){
624  int res=0;
625  //loop over database rows, pushing particles into mctruth object
626  while(1){
627  res = sqlite3_step(statement);
628  if ( res == SQLITE_ROW ){
629  /*
630  * Memo columns:
631  * [0] shower
632  * [1] particle ID (PDG)
633  * [2] momentum: x component [GeV/c]
634  * [3] momentum: y component [GeV/c]
635  * [4] momentum: z component [GeV/c]
636  * [5] position: x component [cm]
637  * [6] position: z component [cm]
638  * [7] time [ns]
639  * [8] energy [GeV]
640  */
641  shower=sqlite3_column_int(statement,0);
642  if(shower!=lastShower){
643  //each new shower gets its own random time and position offsets
644  showerTime=1e9*(flat()*fSampleTime); //converting from s to ns
645  showerTimex=1e9*(flat()*fSampleTime); //converting from s to ns
646  showerTimez=1e9*(flat()*fSampleTime); //converting from s to ns
647  //and a random offset in both z and x controlled by the fRandomXZShift parameter
648  showerXOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
649  showerZOffset=flat()*fRandomXZShift - (fRandomXZShift/2);
650  }
651  pdg=sqlite3_column_int(statement,1);
652  //get mass for this particle
653  double m = 0.; // in GeV
654  TParticlePDG* pdgp = pdgt->GetParticle(pdg);
655  if (pdgp) m = pdgp->Mass();
656 
657  //Note: position/momentum in db have north=-x and west=+z, rotate so that +z is north and +x is west
658  //get momentum components
659  px=sqlite3_column_double(statement,4);//uboone x=Particlez
660  py=sqlite3_column_double(statement,3);
661  pz=-sqlite3_column_double(statement,2);//uboone z=-Particlex
662  etot=sqlite3_column_double(statement,8);
663 
664  //get/calculate position components
665  int boxnoX=0,boxnoZ=0;
666  x=wrapvarBoxNo(sqlite3_column_double(statement,6)+showerXOffset,fShowerBounds[0],fShowerBounds[1],boxnoX);
667  z=wrapvarBoxNo(-sqlite3_column_double(statement,5)+showerZOffset,fShowerBounds[4],fShowerBounds[5],boxnoZ);
668  tParticleTime=sqlite3_column_double(statement,7); //time offset, includes propagation time from top of atmosphere
669  //actual particle time is particle surface arrival time
670  //+ shower start time
671  //+ global offset (fcl parameter, in s)
672  //- propagation time through atmosphere
673  //+ boxNo{X,Z} time offset to make grid boxes have different shower times
674  t=tParticleTime+showerTime+(1e9*fToffset)-fToffset_corsika + showerTimex*boxnoX + showerTimez*boxnoZ;
675  //wrap surface arrival so that it's in the desired time window
676  t=wrapvar(t,(1e9*fToffset),1e9*(fToffset+fSampleTime));
677 
678  simb::MCParticle p(ntotalCtr,pdg,"primary",-200,m,1);
679 
680  //project back to wordvol/fProjectToHeight
681  /*
682  * This back propagation goes from a point on the upper surface of
683  * the cryostat back to the edge of the world, except that that
684  * world is cut short by `fProjectToHeight` (`y2`) ceiling.
685  * The projection will most often lie on that ceiling, but it may
686  * end up instead on one of the side edges of the world, or even
687  * outside it.
688  */
689  double xyzo[3];
690  double x0[3]={x,fShowerBounds[3],z};
691  double dx[3]={px,py,pz};
692  this->ProjectToBoxEdge(x0, dx, x1, x2, y1, y2, z1, z2, xyzo);
693 
694  TLorentzVector pos(xyzo[0],xyzo[1],xyzo[2],t);// time needs to be in ns to match GENIE, etc
695  TLorentzVector mom(px,py,pz,etot);
696  p.AddTrajectoryPoint(pos,mom);
697  mctruth.Add(p);
698  ntotalCtr++;
699  lastShower=shower;
700  }else if ( res == SQLITE_DONE ){
701  break;
702  }else{
703  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
704  }
705  }
706  }else{
707  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kthisStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
708  }
709  nShowerCntr=nShowerCntr-nShowerQry;
710  }
711  }
712  }
process_name opflash particleana ie ie ie z
double fToffset
Time offset of sample, defaults to zero (no offset) [s].
double wrapvar(const double var, const double low, const double high)
double fSampleTime
Duration of sample [s].
var pdg
Definition: selectors.fcl:14
double fProjectToHeight
Height to which particles will be projected [cm].
process_name opflash particleana ie x
BEGIN_PROLOG pz
Definition: filemuons.fcl:10
pdgs p
Definition: selectors.fcl:22
process_name shower
Definition: cheaterreco.fcl:51
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
int fShowerInputs
Number of shower inputs to process from.
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
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[])
Propagates a point back to the surface of a box.
CLHEP::HepRandomEngine & fPoisEngine
double wrapvarBoxNo(const double var, const double low, const double high, int &boxno)
std::vector< int > fMaxShowers
double fToffset_corsika
Timing offset to account for propagation time through atmosphere, populated from db.
BEGIN_PROLOG px
Definition: filemuons.fcl:10
process_name physics producers generator physics producers generator physics producers generator py
double fRandomXZShift
Each shower will be shifted by a random amount in xz so that showers won&#39;t repeatedly sample the same...
double fShowerBounds[6]
Boundaries of area over which showers are to be distributed (x(min), x(max), unused, y, z(min), z(max) )
std::vector< double > fNShowersPerEvent
Number of showers to put in each event of duration fSampleTime; one per showerinput.
CLHEP::HepRandomEngine & fGenEngine
void evgen::CORSIKAGen::openDBs ( std::string const &  module_label)
private

Definition at line 372 of file CORSIKAGen_module.cc.

372  {
373  //choose files based on fShowerInputFiles, copy them with ifdh, open them
374  // for c2: statement is unused
375  //sqlite3_stmt *statement;
376  CLHEP::RandFlat flat(fGenEngine);
377 
378  //setup ifdh object
379  if ( ! fIFDH ) fIFDH = new ifdh_ns::ifdh;
380  const char* ifdh_debug_env = std::getenv("IFDH_DEBUG_LEVEL");
381  if ( ifdh_debug_env ) {
382  mf::LogInfo("CORSIKAGen") << "IFDH_DEBUG_LEVEL: " << ifdh_debug_env<<"\n";
383  fIFDH->set_debug(ifdh_debug_env);
384  }
385 
386  //get ifdh path for each file in fShowerInputFiles, put into selectedflist
387  //if 1 file returned, use that file
388  //if >1 file returned, randomly select one file
389  //if 0 returned, make exeption for missing files
390  std::vector<std::pair<std::string,long>> selectedflist;
391  for(int i=0; i<fShowerInputs; i++){
392  if(fShowerInputFiles[i].find("*")==std::string::npos){
393  //if there are no wildcards, don't call findMatchingFiles
394  selectedflist.push_back(std::make_pair(fShowerInputFiles[i],0));
395  mf::LogInfo("CorsikaGen") << "Selected"<<selectedflist.back().first<<"\n";
396  }else{
397  //use findMatchingFiles
398  //default to IFDH approach when wildcards in path
399  std::vector<std::pair<std::string,long>> flist;
400  std::string path(gSystem->DirName(fShowerInputFiles[i].c_str()));
401  std::string pattern(gSystem->BaseName(fShowerInputFiles[i].c_str()));
402  flist = fIFDH->findMatchingFiles(path,pattern);
403  unsigned int selIndex=-1;
404  if(flist.size()==1){ //0th element is the search path:pattern
405  selIndex=0;
406  }else if(flist.size()>1){
407  selIndex= (unsigned int) (flat()*(flist.size()-1)+0.5); //rnd with rounding, dont allow picking the 0th element
408  }else{
409  throw cet::exception("CORSIKAGen") << "No files returned for path:pattern: "<<path<<":"<<pattern<<std::endl;
410  }
411  selectedflist.push_back(flist[selIndex]);
412  mf::LogInfo("CorsikaGen") << "For "<<fShowerInputFiles[i]<<":"<<pattern
413  <<"\nFound "<< flist.size() << " candidate files"
414  <<"\nChoosing file number "<< selIndex << "\n"
415  <<"\nSelected "<<selectedflist.back().first<<"\n";
416  }
417 
418  }
419 
420  //do the fetching, store local filepaths in locallist
421  std::vector<std::string> locallist;
422  for(unsigned int i=0; i<selectedflist.size(); i++){
423  mf::LogInfo("CorsikaGen")
424  << "Fetching: "<<selectedflist[i].first<<" "<<selectedflist[i].second<<"\n";
425  if (fShowerCopyType == "IFDH") {
426  std::string fetchedfile(fIFDH->fetchInput(selectedflist[i].first));
427  MF_LOG_DEBUG("CorsikaGen") << "Fetched; local path: "<<fetchedfile << "; copy type: IFDH";
428  locallist.push_back(fetchedfile);
429  }
430  else if (fShowerCopyType == "DIRECT") {
431  std::string fetchedfile(selectedflist[i].first);
432  MF_LOG_DEBUG("CorsikaGen") << "Fetched; local path: "<<fetchedfile << "; copy type: DIRECT";
433  locallist.push_back(fetchedfile);
434  }
435  else throw cet::exception("CORSIKAGen") << "Error copying shower db: ShowerCopyType must be \"IFDH\" or \"DIRECT\"\n";
436  }
437 
438  //open the files in fShowerInputFilesLocalPaths with sqlite3
439  for(unsigned int i=0; i<locallist.size(); i++){
440  //prepare and execute statement to attach db file
441  int res=sqlite3_open(locallist[i].c_str(),&fdb[i]);
442  if (res!= SQLITE_OK)
443  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";
444  else
445  mf::LogInfo("CORSIKAGen")<<"Attached db "<< locallist[i]<<"\n";
446  }
447  }
ifdh_ns::ifdh * fIFDH
(optional) flux file handling
std::vector< std::string > fShowerInputFiles
Set of CORSIKA shower data files to use.
BEGIN_PROLOG triggeremu_data_config_icarus settings PMTADCthresholds sequence::icarus_stage0_multiTPC_TPC physics sequence::icarus_stage0_EastHits_TPC physics sequence::icarus_stage0_WestHits_TPC physics producers purityana0 caloskimCalorimetryCryoE physics caloskimCalorimetryCryoW physics path
int fShowerInputs
Number of shower inputs to process from.
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
std::string fShowerCopyType
CLHEP::HepRandomEngine & fGenEngine
void evgen::CORSIKAGen::populateNShowers ( )
private

Definition at line 487 of file CORSIKAGen_module.cc.

487  {
488  //populate vector of the number of showers per event based on:
489  //AREA the showers are being distributed over
490  //TIME of the event (fSampleTime)
491  //flux constants that determine the overall normalizations (fShowerFluxConstants)
492  //Energy range over which the sample was generated (ERANGE_*)
493  //power spectrum over which the sample was generated (ESLOPE)
494 
495 
496  //compute shower area based on the maximal x,z dimensions of cryostat boundaries + fShowerAreaExtension
497  art::ServiceHandle<geo::Geometry const> geom;
498  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
499  double bounds[6] = {0.};
500  geom->CryostatBoundaries(bounds, c);
501  for (unsigned int bnd = 0; bnd<6; bnd++){
502  mf::LogVerbatim("CORSIKAGen")<<"Cryo Boundary: "<<bnd<<"="<<bounds[bnd]<<" ( + Buffer="<<fBuffBox[bnd]<<")\n";
503  if(fabs(bounds[bnd])>fabs(fShowerBounds[bnd])){
504  fShowerBounds[bnd]=bounds[bnd];
505  }
506  }
507  }
508  //add on fShowerAreaExtension without being clever
513 
514  double showersArea=(fShowerBounds[1]/100-fShowerBounds[0]/100)*(fShowerBounds[5]/100-fShowerBounds[4]/100);
515 
516  mf::LogInfo("CORSIKAGen")
517  << "Area extended by : "<<fShowerAreaExtension
518  <<"\nShowers to be distributed betweeen: x="<<fShowerBounds[0]<<","<<fShowerBounds[1]
519  <<" & z="<<fShowerBounds[4]<<","<<fShowerBounds[5]
520  <<"\nShowers to be distributed with random XZ shift: "<<fRandomXZShift<<" cm"
521  <<"\nShowers to be distributed over area: "<<showersArea<<" m^2"
522  <<"\nShowers to be distributed over time: "<<fSampleTime<<" s"
523  <<"\nShowers to be distributed with time offset: "<<fToffset<<" s"
524  <<"\nShowers to be distributed at y: "<<fShowerBounds[3]<<" cm"
525  ;
526 
527  //db variables
528  sqlite3_stmt *statement;
529  const std::string kStatement("select erange_high,erange_low,eslope,nshow from input");
530  double upperLimitOfEnergyRange=0.,lowerLimitOfEnergyRange=0.,energySlope=0.,oneMinusGamma=0.,EiToOneMinusGamma=0.,EfToOneMinusGamma=0.;
531 
532  for(int i=0; i<fShowerInputs; i++){
533  //build and do query to get run info from databases
534  // double thisrnd=flat();//need a new random number for each query
535  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
536  int res=0;
537  res = sqlite3_step(statement);
538  if ( res == SQLITE_ROW ){
539  upperLimitOfEnergyRange=sqlite3_column_double(statement,0);
540  lowerLimitOfEnergyRange=sqlite3_column_double(statement,1);
541  energySlope = sqlite3_column_double(statement,2);
542  fMaxShowers.push_back(sqlite3_column_int(statement,3));
543  oneMinusGamma = 1 + energySlope;
544  EiToOneMinusGamma = pow(lowerLimitOfEnergyRange, oneMinusGamma);
545  EfToOneMinusGamma = pow(upperLimitOfEnergyRange, oneMinusGamma);
546  mf::LogVerbatim("CORSIKAGen")<<"For showers input "<< i<<" found e_hi="<<upperLimitOfEnergyRange<<", e_lo="<<lowerLimitOfEnergyRange<<", slope="<<energySlope<<", k="<<fShowerFluxConstants[i]<<"\n";
547  }else{
548  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
549  }
550  }else{
551  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
552  }
553 
554  //this is computed, how?
555  double NShowers=( M_PI * showersArea * fShowerFluxConstants[i] * (EfToOneMinusGamma - EiToOneMinusGamma) / oneMinusGamma )*fSampleTime;
556  fNShowersPerEvent.push_back(NShowers);
557  mf::LogVerbatim("CORSIKAGen")<<"For showers input "<< i
558  <<" the number of showers per event is "<<(int)NShowers<<"\n";
559  }
560  }
double fToffset
Time offset of sample, defaults to zero (no offset) [s].
double fSampleTime
Duration of sample [s].
int fShowerInputs
Number of shower inputs to process from.
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
std::vector< double > fShowerFluxConstants
Set of flux constants to be associated with each shower data file.
std::vector< int > fMaxShowers
double fRandomXZShift
Each shower will be shifted by a random amount in xz so that showers won&#39;t repeatedly sample the same...
double fShowerBounds[6]
Boundaries of area over which showers are to be distributed (x(min), x(max), unused, y, z(min), z(max) )
std::vector< double > fNShowersPerEvent
Number of showers to put in each event of duration fSampleTime; one per showerinput.
double fShowerAreaExtension
Extend distribution of corsika particles in x,z by this much (e.g. 1000 will extend 10 m in -x...
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 evgen::CORSIKAGen::populateTOffset ( )
private

Definition at line 460 of file CORSIKAGen_module.cc.

460  {
461  //populate TOffset_corsika by finding minimum ParticleTime from db file
462 
463  sqlite3_stmt *statement;
464  const std::string kStatement("select min(t) from particles");
465  double t=0.;
466 
467  for(int i=0; i<fShowerInputs; i++){
468  //build and do query to get run min(t) from each db
469  if ( sqlite3_prepare(fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
470  int res=0;
471  res = sqlite3_step(statement);
472  if ( res == SQLITE_ROW ){
473  t=sqlite3_column_double(statement,0);
474  mf::LogInfo("CORSIKAGen")<<"For showers input "<< i<<" found particles.min(t)="<<t<<"\n";
475  if (i==0 || t<fToffset_corsika) fToffset_corsika=t;
476  }else{
477  throw cet::exception("CORSIKAGen") << "Unexpected sqlite3_step return value: (" <<res<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
478  }
479  }else{
480  throw cet::exception("CORSIKAGen") << "Error preparing statement: (" <<kStatement<<"); "<<"ERROR:"<<sqlite3_errmsg(fdb[i])<<"\n";
481  }
482  }
483 
484  mf::LogInfo("CORSIKAGen")<<"Found corsika timeoffset [ns]: "<< fToffset_corsika<<"\n";
485  }
int fShowerInputs
Number of shower inputs to process from.
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
double fToffset_corsika
Timing offset to account for propagation time through atmosphere, populated from db.
void evgen::CORSIKAGen::produce ( art::Event &  evt)

Definition at line 720 of file CORSIKAGen_module.cc.

720  {
721  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
722 
723  art::ServiceHandle<geo::Geometry const> geom;
724 
725  simb::MCTruth truth;
726  truth.SetOrigin(simb::kCosmicRay);
727 
728  simb::MCTruth pretruth;
729  GetSample(pretruth);
730  mf::LogInfo("CORSIKAGen")<<"GetSample number of particles returned: "<<pretruth.NParticles()<<"\n";
731  // loop over particles in the truth object
732  for(int i = 0; i < pretruth.NParticles(); ++i){
733  simb::MCParticle particle = pretruth.GetParticle(i);
734  const TLorentzVector& v4 = particle.Position();
735  const TLorentzVector& p4 = particle.Momentum();
736  double x0[3] = {v4.X(), v4.Y(), v4.Z() };
737  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
738 
739  // now check if the particle goes through any cryostat in the detector
740  // if so, add it to the truth object.
741  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
742  double bounds[6] = {0.};
743  geom->CryostatBoundaries(bounds, c);
744 
745  //add a buffer box around the cryostat bounds to increase the acceptance and account for scattering
746  //By default, the buffer box has zero size
747  for (unsigned int cb=0; cb<6; cb++)
748  bounds[cb] = bounds[cb]+fBuffBox[cb];
749 
750  //calculate the intersection point with each cryostat surface
751  bool intersects_cryo = false;
752  for (int bnd=0; bnd!=6; ++bnd) {
753  if (bnd<2) {
754  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])};
755  if ( p2[1] >= bounds[2] && p2[1] <= bounds[3] &&
756  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
757  intersects_cryo = true;
758  break;
759  }
760  }
761  else if (bnd>=2 && bnd<4) {
762  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])};
763  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
764  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
765  intersects_cryo = true;
766  break;
767  }
768  }
769  else if (bnd>=4) {
770  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]};
771  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
772  p2[1] >= bounds[2] && p2[1] <= bounds[3] ) {
773  intersects_cryo = true;
774  break;
775  }
776  }
777  }
778 
779  if (intersects_cryo){
780  truth.Add(particle);
781  break; //leave loop over cryostats to avoid adding particle multiple times
782  }// end if particle goes into a cryostat
783  }// end loop over cryostats in the detector
784 
785  }// loop on particles
786 
787  mf::LogInfo("CORSIKAGen")<<"Number of particles from getsample crossing cryostat + bounding box: "<<truth.NParticles()<<"\n";
788 
789  truthcol->push_back(truth);
790  evt.put(std::move(truthcol));
791 
792  return;
793  }// end produce
const Cut kCosmicRay
TCEvent evt
Definition: DataStructs.cxx:8
void GetSample(simb::MCTruth &)
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 evgen::CORSIKAGen::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[] 
)
private

Propagates a point back to the surface of a box.

Parameters
xyzcoordinates of the point to be propagated ({ x, y, z })
dxyzdirection of the point ({ dx, dy, dz })
xlolower x coordinate of the target box
xhiupper x coordinate of the target box
ylolower y coordinate of the target box
yhiupper y coordinate of the target box
zlolower z coordinate of the target box
zhiupper z coordinate of the target box
xyzout_(output, room for at least 3 numbers)_ propagated point

The point xyz, assumed to be inside the box, is propagated at the level of the closest among the sides of the box. Note that this means the propagated point might still be not on the surface of the box, even if it would later reach it. It is anyway guaranteed that xyzout is not inside the target box, and that at least one of its three coordinates is on the box surface.

Definition at line 333 of file CORSIKAGen_module.cc.

341  {
342 
343 
344  //we want to project backwards, so take mirror of momentum
345  const double dxyz[3]={-indxyz[0],-indxyz[1],-indxyz[2]};
346 
347  // Compute the distances to the x/y/z walls
348  double dx = 99.E99;
349  double dy = 99.E99;
350  double dz = 99.E99;
351  if (dxyz[0] > 0.0) { dx = (xhi-xyz[0])/dxyz[0]; }
352  else if (dxyz[0] < 0.0) { dx = (xlo-xyz[0])/dxyz[0]; }
353  if (dxyz[1] > 0.0) { dy = (yhi-xyz[1])/dxyz[1]; }
354  else if (dxyz[1] < 0.0) { dy = (ylo-xyz[1])/dxyz[1]; }
355  if (dxyz[2] > 0.0) { dz = (zhi-xyz[2])/dxyz[2]; }
356  else if (dxyz[2] < 0.0) { dz = (zlo-xyz[2])/dxyz[2]; }
357 
358 
359  // Choose the shortest distance
360  double d = 0.0;
361  if (dx < dy && dx < dz) d = dx;
362  else if (dy < dz && dy < dx) d = dy;
363  else if (dz < dx && dz < dy) d = dz;
364 
365  // Make the step
366  for (int i = 0; i < 3; ++i) {
367  xyzout[i] = xyz[i] + dxyz[i]*d;
368  }
369 
370  }
double evgen::CORSIKAGen::wrapvar ( const double  var,
const double  low,
const double  high 
)
private

Definition at line 449 of file CORSIKAGen_module.cc.

449  {
450  //wrap variable so that it's always between low and high
451  return (var - (high - low) * floor(var/(high-low))) + low;
452  }
standard_dbscan3dalg useful for diagnostics hits not in a line will not be clustered on on only for track like only for track like on on the smaller the less shower like tracks low
double evgen::CORSIKAGen::wrapvarBoxNo ( const double  var,
const double  low,
const double  high,
int &  boxno 
)
private

Definition at line 454 of file CORSIKAGen_module.cc.

454  {
455  //wrap variable so that it's always between low and high
456  boxno=int(floor(var/(high-low)));
457  return (var - (high - low) * floor(var/(high-low))) + low;
458  }
standard_dbscan3dalg useful for diagnostics hits not in a line will not be clustered on on only for track like only for track like on on the smaller the less shower like tracks low

Member Data Documentation

std::vector<double> evgen::CORSIKAGen::fBuffBox
private

Buffer box extensions to cryostat in each direction (6 of them: x_lo,x_hi,y_lo,y_hi,z_lo,z_hi) [cm].

Definition at line 280 of file CORSIKAGen_module.cc.

sqlite3* evgen::CORSIKAGen::fdb[5]
private

Pointers to sqlite3 database object, max of 5.

Definition at line 282 of file CORSIKAGen_module.cc.

CLHEP::HepRandomEngine& evgen::CORSIKAGen::fGenEngine
private

Definition at line 284 of file CORSIKAGen_module.cc.

ifdh_ns::ifdh* evgen::CORSIKAGen::fIFDH =0
private

(optional) flux file handling

Definition at line 271 of file CORSIKAGen_module.cc.

std::vector<int> evgen::CORSIKAGen::fMaxShowers
private

Definition at line 268 of file CORSIKAGen_module.cc.

std::vector<double> evgen::CORSIKAGen::fNShowersPerEvent
private

Number of showers to put in each event of duration fSampleTime; one per showerinput.

Definition at line 267 of file CORSIKAGen_module.cc.

CLHEP::HepRandomEngine& evgen::CORSIKAGen::fPoisEngine
private

Definition at line 285 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fProjectToHeight =0.
private

Height to which particles will be projected [cm].

Definition at line 274 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fRandomXZShift =0.
private

Each shower will be shifted by a random amount in xz so that showers won't repeatedly sample the same space [cm].

Definition at line 283 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fSampleTime =0.
private

Duration of sample [s].

Definition at line 278 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fShowerAreaExtension =0.
private

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].

Definition at line 281 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fShowerBounds[6] ={0.,0.,0.,0.,0.,0.}
private

Boundaries of area over which showers are to be distributed (x(min), x(max), unused, y, z(min), z(max) )

Definition at line 269 of file CORSIKAGen_module.cc.

std::string evgen::CORSIKAGen::fShowerCopyType
private

Definition at line 276 of file CORSIKAGen_module.cc.

std::vector< double > evgen::CORSIKAGen::fShowerFluxConstants
private

Set of flux constants to be associated with each shower data file.

Definition at line 277 of file CORSIKAGen_module.cc.

std::vector< std::string > evgen::CORSIKAGen::fShowerInputFiles
private

Set of CORSIKA shower data files to use.

Definition at line 275 of file CORSIKAGen_module.cc.

int evgen::CORSIKAGen::fShowerInputs =0
private

Number of shower inputs to process from.

Definition at line 266 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fToffset =0.
private

Time offset of sample, defaults to zero (no offset) [s].

Definition at line 279 of file CORSIKAGen_module.cc.

double evgen::CORSIKAGen::fToffset_corsika =0.
private

Timing offset to account for propagation time through atmosphere, populated from db.

Definition at line 270 of file CORSIKAGen_module.cc.


The documentation for this class was generated from the following file: