1 #ifndef EVGEN_CORSIKAGen_H 8 #define EVGEN_CORSIKAGen_H 12 #include "TDatabasePDG.h" 42 #include "CLHEP/Random/RandFlat.h" 43 #include "CLHEP/Random/RandPoissonQ.h" 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);
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.))
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();
114 if(fSampleTime==0.)
throw cet::exception(
"CORSIKAGen") <<
"SampleTime not set!";
116 if(fProjectToHeight==0.)
mf::LogInfo(
"CORSIKAGen")<<
"Using 0. for fProjectToHeight!" 123 this->reconfigure(p);
126 this->populateNShowers();
127 this->populateTOffset();
129 produces< std::vector<simb::MCTruth> >();
130 produces< sumdata::RunData, art::InRun >();
136 sqlite3_close(
fdb[i]);
143 const double indxyz[],
154 const double dxyz[3]={-indxyz[0],-indxyz[1],-indxyz[2]};
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]; }
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;
175 for (
int i = 0; i < 3; ++i) {
176 xyzout[i] = xyz[i] + dxyz[i]*
d;
187 CLHEP::HepRandomEngine &engine = rng->
getEngine(
"gen");
188 CLHEP::RandFlat flat(engine);
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);
202 std::vector<std::pair<std::string,long>> selectedflist;
207 mf::LogInfo(
"CorsikaGen") <<
"Selected"<<selectedflist.back().first<<
"\n";
210 std::vector<std::pair<std::string,long>> flist;
213 flist =
fIFDH->findMatchingFiles(path,pattern);
214 unsigned int selIndex=-1;
217 }
else if(flist.size()>1){
218 selIndex= (
unsigned int) (flat()*(flist.size()-1)+0.5);
220 throw cet::exception(
"CORSIKAGen") <<
"No files returned for path:pattern: "<<path<<
":"<<pattern<<std::endl;
222 selectedflist.push_back(flist[selIndex]);
224 <<
"\nFound "<< flist.size() <<
" candidate files" 225 <<
"\nChoosing file number "<< selIndex <<
"\n" 226 <<
"\nSelected "<<selectedflist.back().first<<
"\n";
232 std::vector<std::string> locallist;
233 for(
unsigned int i=0; i<selectedflist.size(); i++){
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);
242 for(
unsigned int i=0; i<locallist.size(); i++){
244 int res=sqlite3_open(locallist[i].c_str(),&
fdb[i]);
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";
248 mf::LogInfo(
"CORSIKAGen")<<
"Attached db "<< locallist[i]<<
"\n";
254 return (var - (high - low) * floor(var/(high-low))) + low;
259 boxno=int(floor(var/(high-low)));
260 return (var - (high - low) * floor(var/(high-low))) + low;
266 sqlite3_stmt *statement;
267 const std::string kStatement(
"select min(t) from particles");
272 if ( sqlite3_prepare(
fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
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";
280 throw cet::exception(
"CORSIKAGen") <<
"Unexpected sqlite3_step return value: (" <<res<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
283 throw cet::exception(
"CORSIKAGen") <<
"Error preparing statement: (" <<kStatement<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
301 for(
unsigned int c = 0; c < geom->
Ncryostats(); ++c){
304 for (
unsigned int bnd = 0; bnd<6; bnd++){
320 <<
"Area extended by : "<<fShowerAreaExtension
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" 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.;
345 if ( sqlite3_prepare(
fdb[i], kStatement.c_str(), -1, &statement, 0 ) == SQLITE_OK ){
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";
358 throw cet::exception(
"CORSIKAGen") <<
"Unexpected sqlite3_step return value: (" <<res<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
361 throw cet::exception(
"CORSIKAGen") <<
"Error preparing statement: (" <<kStatement<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
368 <<
" the number of showers per event is "<<(int)NShowers<<
"\n";
385 static TDatabasePDG* pdgt = TDatabasePDG::Instance();
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)");
393 CLHEP::HepRandomEngine &engine = rng->
getEngine(
"gen");
394 CLHEP::RandFlat flat(engine);
396 CLHEP::HepRandomEngine &engine_pois = rng->
getEngine(
"pois");
397 CLHEP::RandPoissonQ randpois(engine_pois);
404 geom->
WorldBox(&x1, &x2, &y1, &y2, &z1, &z2);
408 double fBoxDelta=1.e-5;
422 double px,py,pz,
x,
z,tParticleTime,etot,showerTime=0.,showerTimex=0.,showerTimez=0.,showerXOffset=0.,showerZOffset=0.,t;
427 while(nShowerCntr>0){
432 nShowerQry=nShowerCntr;
435 double thisrnd=flat();
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 ){
442 res = sqlite3_step(statement);
443 if ( res == SQLITE_ROW ){
444 shower=sqlite3_column_int(statement,0);
445 if(shower!=lastShower){
454 pdg=sqlite3_column_int(statement,1);
457 TParticlePDG* pdgp = pdgt->GetParticle(pdg);
458 if (pdgp) m = pdgp->Mass();
462 px=sqlite3_column_double(statement,4);
463 py=sqlite3_column_double(statement,3);
464 pz=-sqlite3_column_double(statement,2);
465 etot=sqlite3_column_double(statement,8);
468 int boxnoX=0,boxnoZ=0;
471 tParticleTime=sqlite3_column_double(statement,7);
486 double dx[3]={px,py,pz};
489 TLorentzVector pos(xyzo[0],xyzo[1],xyzo[2],t);
490 TLorentzVector mom(px,py,pz,etot);
495 }
else if ( res == SQLITE_DONE ){
498 throw cet::exception(
"CORSIKAGen") <<
"Unexpected sqlite3_step return value: (" <<res<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
502 throw cet::exception(
"CORSIKAGen") <<
"Error preparing statement: (" <<kthisStatement<<
"); "<<
"ERROR:"<<sqlite3_errmsg(
fdb[i])<<
"\n";
504 nShowerCntr=nShowerCntr-nShowerQry;
527 run.
put(std::move(runcol));
533 std::unique_ptr< std::vector<simb::MCTruth> > truthcol(
new std::vector<simb::MCTruth>);
544 for(
int i = 0; i < pretruth.
NParticles(); ++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()};
553 for(
unsigned int c = 0; c < geom->
Ncryostats(); ++c){
559 for (
unsigned int cb=0; cb<6; cb++)
560 bounds[cb] = bounds[cb]+
fBuffBox[cb];
563 bool intersects_cryo =
false;
564 for (
int bnd=0; bnd!=6; ++bnd) {
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;
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;
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;
591 if (intersects_cryo){
599 mf::LogInfo(
"CORSIKAGen")<<
"Number of particles from getsample crossing cryostat + bounding box: "<<truth.
NParticles()<<
"\n";
601 truthcol->push_back(truth);
602 evt.
put(std::move(truthcol));
616 #endif // EVGEN_CORSIKAGen_H 617
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
double fSampleTime
Duration of sample [s].
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
double fProjectToHeight
Height to which particles will be projected [cm].
Encapsulate the construction of a single cyostat.
Float_t y1[n_points_granero]
Collect all the geometry header files together.
void SetOrigin(simb::Origin_t origin)
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]
A module to check the results from the Monte Carlo generator.
ifdh_ns::ifdh * fIFDH
(optional) flux file handling
art::ProductID put(std::unique_ptr< PROD > &&)
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
void Add(simb::MCParticle &part)
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.
Float_t y2[n_points_geant4]
int fShowerInputs
Number of shower inputs to process from.
ProductID put(std::unique_ptr< PROD > &&product)
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
base_engine_t & getEngine() const
sqlite3 * fdb[5]
Pointers to sqlite3 database object, max of 5.
#define DEFINE_ART_MODULE(klass)
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[])
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
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
double fRandomXZShift
Each shower will be shifted by a random amount in xz so that showers won't repeatedly sample the same...
double fShowerBounds[6]
Boundaries of area over which showers are to be distributed.
CORSIKAGen(fhicl::ParameterSet const &pset)
void beginRun(art::Run &run)
Event generator information.
Float_t x2[n_points_geant4]
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)
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
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].