LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
MUSUN_module.cc
Go to the documentation of this file.
1 
143 // C++ includes.
144 #include <cmath>
145 #include <cstdlib>
146 #include <exception>
147 #include <fstream>
148 #include <iostream>
149 #include <memory>
150 #include <sstream>
151 #include <string>
152 
153 // Framework includes
159 #include "art_root_io/TFileService.h"
160 #include "cetlib_except/exception.h"
161 #include "fhiclcpp/ParameterSet.h"
163 
164 // art extensions
166 
167 // nusimdata includes
170 
171 // lar includes
174 
175 #include "TDatabasePDG.h"
176 #include "TTree.h"
177 
178 #include "CLHEP/Random/RandFlat.h"
179 #include "CLHEP/Random/RandGaussQ.h"
180 
181 // for c2: NTPCs is no longer used
182 //const int NTPCs = 300;
183 
184 namespace simb {
185  class MCTruth;
186 }
187 
188 namespace evgen {
189 
191  class MUSUN : public art::EDProducer {
192 
193  public:
194  explicit MUSUN(fhicl::ParameterSet const& pset);
195 
196  // This is called for each event.
197  void produce(art::Event& evt);
198  void beginJob();
199  void beginRun(art::Run& run);
200  void endRun(art::Run& run);
201 
202  private:
203  void SampleOne(unsigned int i, simb::MCTruth& mct, CLHEP::HepRandomEngine& engine);
204 
205  void initialization(double theta1,
206  double theta2,
207  double phi1,
208  double phi2,
209  int figflag,
210  double s_hor,
211  double s_ver1,
212  double s_ver2,
213  double& FI);
214 
215  void sampling(double& E,
216  double& theta,
217  double& phi,
218  double& dep,
219  CLHEP::HepRandomEngine& engine);
220 
221  static const int kGAUS = 1;
222 
223  CLHEP::HepRandomEngine& fEngine;
224 
225  int fPDG;
226  double fChargeRatio;
227 
228  std::string fInputDir;
229  std::string fInputFile1;
230  std::string fInputFile2;
231  std::string fInputFile3;
232 
233  double fCavernAngle;
234  double fRockDensity;
235 
239  double fEmin;
240  double fEmax;
241 
242  double fThetamin;
243  double fThetamax;
244  double fPhimin;
245  double fPhimax;
246 
247  int figflag;
248  double fXmin;
249  double fYmin;
250  double fZmin;
251  double fXmax;
252  double fYmax;
253  double fZmax;
254 
257  double
259  double
261  double
263 
264  double fT0;
265  double fSigmaT;
266  int fTDist;
267 
268  //Define TFS histograms.....
269  /*
270  TH1D* hPDGCode;
271  TH1D* hPositionX;
272  TH1D* hPositionY;
273  TH1D* hPositionZ;
274  TH1D* hTime;
275  TH1D* hMomentumHigh;
276  TH1D* hMomentum;
277  TH1D* hEnergyHigh;
278  TH1D* hEnergy;
279  TH1D* hDepth;
280  TH1D* hDirCosineX;
281  TH1D* hDirCosineY;
282  TH1D* hDirCosineZ;
283  TH1D* hTheta;
284  TH1D* hPhi;
285  */
286  int PdgCode;
287  double Energy, phi, theta, dep, Time;
288  double Momentum, px0, py0, pz0;
289  double x0, y0, z0, cx, cy, cz;
290 
291  //Define some variables....
292  double FI = 0.;
293  double s_hor = 0.;
294  double s_ver1 = 0.;
295  double s_ver2 = 0.;
296 
297  double spmu[121][62][51];
298  double fnmu[32401];
299  double depth[360][91];
300  double fmu[360][91];
301  // for c2: e1 and e2 are unused
302  //double e1, e2, the1, the2, ph1, ph2;
303  double the1, the2, ph1, ph2;
304  double se = 0.;
305  double st = 0.;
306  double sp = 0.;
307  double sd = 0.;
308 
309  unsigned int NEvents = 0;
310 
311  // TTree
312  TTree* fTree;
313  /*
314  // TTree for CryoPos
315  TTree* fCryos;
316  int NumTPCs;
317  double TPCMinX[NTPCs];
318  double TPCMaxX[NTPCs];
319  double TPCMinY[NTPCs];
320  double TPCMaxY[NTPCs];
321  double TPCMinZ[NTPCs];
322  double TPCMaxZ[NTPCs];
323  double CryoSize[6];
324  double DetHall[6];
325  */
326  };
327 }
328 
329 namespace evgen {
330 
331  //____________________________________________________________________________
332  MUSUN::MUSUN(fhicl::ParameterSet const& pset)
333  : art::EDProducer{pset}
334  // create a default random engine; obtain the random seed from NuRandomService,
335  // unless overridden in configuration with key "Seed"
336  , fEngine(art::ServiceHandle<rndm::NuRandomService>()->registerAndSeedEngine(createEngine(0),
337  pset,
338  "Seed"))
339  , fPDG{pset.get<int>("PDG")}
340  , fChargeRatio{pset.get<double>("ChargeRatio")}
341  , fInputDir{pset.get<std::string>("InputDir")}
342  , fInputFile1{pset.get<std::string>("InputFile1")}
343  , fInputFile2{pset.get<std::string>("InputFile2")}
344  , fInputFile3{pset.get<std::string>("InputFile3")}
345  , fCavernAngle{pset.get<double>("CavernAngle")}
346  , fRockDensity{pset.get<double>("RockDensity")}
347  , fEmin{pset.get<double>("Emin")}
348  , fEmax{pset.get<double>("Emax")}
349  , fThetamin{pset.get<double>("Thetamin")}
350  , fThetamax{pset.get<double>("Thetamax")}
351  , fPhimin{pset.get<double>("Phimin")}
352  , fPhimax{pset.get<double>("Phimax")}
353  , figflag{pset.get<int>("igflag")}
354  , fXmin{pset.get<double>("Xmin")}
355  , fYmin{pset.get<double>("Ymin")}
356  , fZmin{pset.get<double>("Zmin")}
357  , fXmax{pset.get<double>("Xmax")}
358  , fYmax{pset.get<double>("Ymax")}
359  , fZmax{pset.get<double>("Zmax")}
360  , fDetRotX{pset.get<double>("DetRotX")}
361  , fDetRotY{pset.get<double>("DetRotY")}
362  , fDetRotZ{pset.get<double>("DetRotZ")}
363  , fT0{pset.get<double>("T0")}
364  , fSigmaT{pset.get<double>("SigmaT")}
365  , fTDist{pset.get<int>("TDist")}
366  {
367  produces<std::vector<simb::MCTruth>>();
368  produces<sumdata::RunData, art::InRun>();
369  }
370 
372  // Begin Job
375  {
376 
377  fDetRotX = fDetRotX * M_PI / 180.;
378  fDetRotY = fDetRotY * M_PI / 180.;
379  fDetRotZ = fDetRotZ * M_PI / 180.;
380 
381  MF_LOG_WARNING("MUSUN") << "Rotation: X = " << fDetRotX << ", Y = " << fDetRotY
382  << ", Z = " << fDetRotZ << std::endl;
383 
384  // Make the Histograms....
386  /*
387  hPDGCode = tfs->make<TH1D>("hPDGCode" ,"PDG Code" ,30 , -15 , 15 );
388  hPositionX = tfs->make<TH1D>("hPositionX" ,"Position (cm)" ,500, ( fXmin - 10 ), ( fXmax + 10 ) );
389  hPositionY = tfs->make<TH1D>("hPositionY" ,"Position (cm)" ,500, ( fYmin - 10 ), ( fYmax + 10 ) );
390  hPositionZ = tfs->make<TH1D>("hPositionZ" ,"Position (cm)" ,500, ( fZmin - 10 ), ( fZmax + 10 ) );
391  hTime = tfs->make<TH1D>("hTime" ,"Time (s)" ,500, 0 , 1e6 );
392  hMomentumHigh = tfs->make<TH1D>("hMomentumHigh","Momentum (GeV)",500, fEmin, fEmax );
393  hMomentum = tfs->make<TH1D>("hMomentum" ,"Momentum (GeV)",500, fEmin, fEmin+1e3 );
394  hEnergyHigh = tfs->make<TH1D>("hEnergyHigh" ,"Energy (GeV)" ,500, fEmin, fEmax );
395  hEnergy = tfs->make<TH1D>("hEnergy" ,"Energy (GeV)" ,500, fEmin, fEmin+1e3 );
396  hDepth = tfs->make<TH1D>("hDepth" ,"Depth (m)" ,800, 0 , 14000 );
397 
398  hDirCosineX = tfs->make<TH1D>("hDirCosineX","Normalised Direction cosine",500, -1, 1 );
399  hDirCosineY = tfs->make<TH1D>("hDirCosineY","Normalised Direction cosine",500, -1, 1 );
400  hDirCosineZ = tfs->make<TH1D>("hDirCosineZ","Normalised Direction cosine",500, -1, 1 );
401 
402  hTheta = tfs->make<TH1D>("hTheta" ,"Angle (degrees)",500, 0, 90 );
403  hPhi = tfs->make<TH1D>("hPhi" ,"Angle (degrees)",500, 0, 365 );
404  */
405  fTree = tfs->make<TTree>("Generator", "analysis tree");
406  fTree->Branch("particleID", &PdgCode, "particleID/I");
407  fTree->Branch("energy", &Energy, "energy/D");
408  fTree->Branch("time", &Time, "Time/D");
409  fTree->Branch("posX", &x0, "posX/D");
410  fTree->Branch("posY", &y0, "posY/D");
411  fTree->Branch("posZ", &z0, "posZ/D");
412  fTree->Branch("cosX", &cx, "cosX/D");
413  fTree->Branch("cosY", &cy, "cosY/D");
414  fTree->Branch("cosZ", &cz, "cosZ/D");
415  fTree->Branch("theta", &theta, "theta/D");
416  fTree->Branch("phi", &phi, "phi/D");
417  fTree->Branch("depth", &dep, "dep/D");
418  /*
419  fCryos = tfs->make<TTree>("CryoSizes","cryo tree");
420  fCryos->Branch("NumTPCs" , &NumTPCs , "NumTPCs/I" );
421  fCryos->Branch("TPCMinX" , &TPCMinX , "TPCMinX[NumTPCs]/D");
422  fCryos->Branch("TPCMaxX" , &TPCMaxX , "TPCMaxX[NumTPCs]/D");
423  fCryos->Branch("TPCMinY" , &TPCMinY , "TPCMinY[NumTPCs]/D");
424  fCryos->Branch("TPCMaxY" , &TPCMaxY , "TPCMaxY[NumTPCs]/D");
425  fCryos->Branch("TPCMinZ" , &TPCMinZ , "TPCMinZ[NumTPCs]/D");
426  fCryos->Branch("TPCMaxZ" , &TPCMaxZ , "TPCMaxZ[NumTPCs]/D");
427  fCryos->Branch("CryoSize", &CryoSize, "CryoSize[6]/D" );
428  fCryos->Branch("DetHall" , &DetHall , "DetHall[6]/D" );
429  */
430  }
431 
433  // Begin Run
436  {
437  // Check fcl parameters were set correctly
438  if (fThetamax > 90.5)
439  throw cet::exception("MUSUNGen")
440  << "\nThetamax has to be less than " << M_PI / 2 << ", but was entered as " << fThetamax
441  << ", this causes an error so leaving program now...\n\n";
442  if (fThetamin < 0)
443  throw cet::exception("MUSUNGen")
444  << "\nThetamin has to be more than 0, but was entered as " << fThetamin
445  << ", this causes an error so leaving program now...\n\n";
446  if (fThetamax < fThetamin)
447  throw cet::exception("MUSUNGen") << "\nMinimum angle is bigger than maximum angle....causes "
448  "an error so leaving program now....\n\n";
449  if (fPhimax > 360.5)
450  throw cet::exception("MUSUNGen")
451  << "\nPhimax has to be less than " << 2 * M_PI << ", but was entered as " << fPhimax
452  << ", this cause an error so leaving program now...\n\n";
453  if (fPhimin < 0)
454  throw cet::exception("MUSUNGen")
455  << "\nPhimin has to be more than 0, but was entered as " << fPhimin
456  << ", this causes an error so leaving program now...\n\n";
457  if (fPhimax < fPhimin)
458  throw cet::exception("MUSUNGen") << "\nMinimum angle is bigger than maximum angle....causes "
459  "an error so leaving program now....\n\n";
460  if (fEmax < fEmin)
461  throw cet::exception("MUSUNGen")
462  << "\nMinimum energy is bigger than maximum energy....causes an error so leaving program "
463  "now....\n\n";
464 
465  // grab the geometry object to see what geometry we are using
467  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()), art::fullRun());
468 
469  // area of the horizontal plane of the parallelepiped
470  s_hor = (fZmax - fZmin) * (fXmax - fXmin);
471  // area of the vertical plane of the parallelepiped, perpendicular to z-axis
472  s_ver1 = (fXmax - fXmin) * (fYmax - fYmin);
473  // area of the vertical plane of the parallelepiped, perpendicular to x-axis
474  s_ver2 = (fZmax - fZmin) * (fYmax - fYmin);
475 
476  //std::cout << s_hor << " " << s_ver1 << " " << s_ver2 << std::endl;
477 
479 
480  std::cout << "Material - SURF rock" << std::endl;
481  std::cout << "Density = " << fRockDensity << " g/cm^3" << std::endl;
482  std::cout << "Parameters for muon spectrum are from LVD best fit" << std::endl;
483  std::cout << "Muon energy range = " << fEmin << " - " << fEmax << " GeV" << std::endl;
484  std::cout << "Zenith angle range = " << fThetamin << " - " << fThetamax << " degrees"
485  << std::endl;
486  std::cout << "Azimuthal angle range = " << fPhimin << " - " << fPhimax << " degrees"
487  << std::endl;
488  std::cout << "Global intensity = " << FI
489  << " (cm^2 s)^(-1) or s^(-1) (for muons on the surface)" << std::endl;
490  /*
491  NumTPCs = geo->NTPC(0);
492  std::cout << "There are " << NumTPCs << " in cryostat 0" << std::endl;
493  for (unsigned int c=0; c<geo->Ncryostats(); c++) {
494  const geo::CryostatGeo& cryostat=geo->Cryostat(c);
495  geo->CryostatBoundaries( CryoSize, 0 );
496  std::cout << "Cryo bounds " << CryoSize[0] << " "<< CryoSize[1] << " "<< CryoSize[2] << " "<< CryoSize[3] << " "<< CryoSize[4] << " "<< CryoSize[5] << std::endl;
497  for (unsigned int t=0; t<cryostat.NTPC(); t++) {
498  geo::TPCID id;
499  id.Cryostat=c;
500  id.TPC=t;
501  id.isValid=true;
502  const geo::TPCGeo& tpc=cryostat.TPC(id);
503  TPCMinX[t] = tpc.MinX();
504  TPCMaxX[t] = tpc.MaxX();
505  TPCMinY[t] = tpc.MinY();
506  TPCMaxY[t] = tpc.MaxY();
507  TPCMinZ[t] = tpc.MinZ();
508  TPCMaxZ[t] = tpc.MaxZ();
509  std::cout << t << "\t" << TPCMinX[t] << " " << TPCMaxX[t] << " " << TPCMinY[t] << " " << TPCMaxY[t] << " " << TPCMinZ[t] << " " << TPCMaxZ[t] << std::endl;
510  }
511  }
512  fCryos -> Fill();
513  */
514  return;
515  }
516 
518  // End Run
521  {
522  std::cout << "\n\nNumber of muons = " << NEvents << std::endl;
523  std::cout << "Mean muon energy = " << se / NEvents << " GeV" << std::endl;
524  std::cout << "Mean zenith angle (theta) = " << st / NEvents << " degrees" << std::endl;
525  std::cout << "Mean azimuthal angle (phi)= " << sp / NEvents << " degrees" << std::endl;
526  std::cout << "Mean slant depth = " << sd / NEvents << " m w.e." << std::endl;
527  }
529  // Produce
532  {
534  std::unique_ptr<std::vector<simb::MCTruth>> truthcol(new std::vector<simb::MCTruth>);
535 
536  simb::MCTruth truth;
538 
539  ++NEvents;
540 
541  SampleOne(NEvents, truth, fEngine);
542 
543  MF_LOG_DEBUG("MUSUN") << truth;
544 
545  truthcol->push_back(truth);
546 
547  evt.put(std::move(truthcol));
548  }
550  // Sample One
552  // Draw the type, momentum and position of a single particle from the
553  // FCIHL description
554  void MUSUN::SampleOne(unsigned int i, simb::MCTruth& mct, CLHEP::HepRandomEngine& engine)
555  {
556  CLHEP::RandFlat flat(engine);
557  CLHEP::RandGaussQ gauss(engine);
558 
559  Energy = 0;
560  theta = 0;
561  phi = 0;
562  dep = 0;
563  Time = 0;
564 
565  sampling(Energy, theta, phi, dep, engine);
566 
567  theta = theta * M_PI / 180;
568 
569  // changing the angle phi so z-axis is positioned along the long side
570  // of the cavern pointing at 14 deg from the North to the East.
571  // phi += (90. - 14.0);
572  // Want our co-ord rotation going from East to South.
573  phi += fCavernAngle;
574  if (phi >= 360.) phi -= 360.;
575  if (phi < 0) phi += 360.;
576  phi *= M_PI / 180.;
577 
578  // set track id to -i as these are all primary particles and have id <= 0
579  int trackid = -1 * (i + 1);
580  std::string primary("primary");
581 
582  // Work out whether particle/antiparticle, and mass...
583  double m = 0.0;
584  PdgCode = fPDG;
585  double ChargeCheck = 1. / (1 + fChargeRatio);
586  double pdgfire = flat.fire();
587  if (pdgfire < ChargeCheck) PdgCode = -PdgCode;
588 
589  static TDatabasePDG pdgt;
590  TParticlePDG* pdgp = pdgt.GetParticle(PdgCode);
591  if (pdgp) m = pdgp->Mass();
592 
593  //std::cout << pdgfire << " " << ChargeCheck << " " << PdgCode << " " << m << std::endl;
594 
595  // Work out T0...
596  if (fTDist == kGAUS) { Time = gauss.fire(fT0, fSigmaT); }
597  else {
598  Time = fT0 + fSigmaT * (2.0 * flat.fire() - 1.0);
599  }
600 
601  // The minus sign above is for y-axis pointing up, so the y-momentum
602  // is always pointing down
603  cx = +sin(theta) * sin(phi);
604  cy = -cos(theta);
605  cz = +sin(theta) * cos(phi);
606  Momentum = std::sqrt(Energy * Energy - m * m); // Get momentum
607  px0 = Momentum * cx;
608  py0 = Momentum * cy;
609  pz0 = Momentum * cz;
610  TLorentzVector pvec(px0, py0, pz0, Energy);
611 
612  // transform into actual detector's coordinates
613  pvec.RotateX(-fDetRotX);
614  pvec.RotateY(-fDetRotY);
615  pvec.RotateZ(-fDetRotZ);
616 
617  // Muon coordinates
618  double sh1 = s_hor * cos(theta);
619  double sv1 = s_ver1 * sin(theta) * fabs(cos(phi));
620  double sv2 = s_ver2 * sin(theta) * fabs(sin(phi));
621  double ss = sh1 + sv1 + sv2;
622  double xfl1 = flat.fire();
623  if (xfl1 <= sh1 / ss) {
624  x0 = (fXmax - fXmin) * flat.fire() + fXmin;
625  y0 = fYmax;
626  z0 = (fZmax - fZmin) * flat.fire() + fZmin;
627  }
628  else if (xfl1 <= (sh1 + sv1) / ss) {
629  x0 = (fXmax - fXmin) * flat.fire() + fXmin;
630  y0 = (fYmax - fYmin) * flat.fire() + fYmin;
631  if (cz >= 0)
632  z0 = fZmin;
633  else
634  z0 = fZmax;
635  }
636  else {
637  if (cx >= 0)
638  x0 = fXmin;
639  else
640  x0 = fXmax;
641  y0 = (fYmax - fYmin) * flat.fire() + fYmin;
642  z0 = (fZmax - fZmin) * flat.fire() + fZmin;
643  }
644  // Make Lorentz vector for x and t....
645  TLorentzVector pos(x0, y0, z0, Time);
646 
647  // transform into actual detector's coordinates
648  pos.RotateX(-fDetRotX);
649  pos.RotateY(-fDetRotY);
650  pos.RotateZ(-fDetRotZ);
651 
652  // Parameters written to the file muons_surf_v2_test*.dat
653  // nmu - muon sequential number
654  // id_part - muon charge (10 - positive, 11 - negative )
655  // Energy - total muon energy in GeV assuming ultrarelativistic muons
656  // x0, y0, z0 - muon coordinates on the surface of parallelepiped
657  // specified above; x-axis and y-axis are pointing in the
658  // directions such that the angle phi (from the slant depth
659  // distribution files) is measured from x to y. z-axis is
660  // pointing upwards.
661  // cx, cy, cz - direction cosines.
662 
663  simb::MCParticle part(trackid, PdgCode, primary);
664  part.AddTrajectoryPoint(pos, pvec);
665 
666  mct.Add(part);
667 
668  theta = theta * 180 / M_PI;
669  phi = phi * 180 / M_PI;
670 
671  // Sum energies, angles, depth for average outputs.
672  se += Energy;
673  st += theta;
674  sp += phi;
675  sd += dep;
676 
677  // Fill Histograms.....
678  /*
679  hPDGCode ->Fill (PdgCode);
680  hPositionX ->Fill (x0);
681  hPositionY ->Fill (y0);
682  hPositionZ ->Fill (z0);
683  hTime ->Fill (Time);
684  hMomentumHigh ->Fill (Momentum);
685  hMomentum ->Fill (Momentum);
686  hEnergyHigh ->Fill (Energy);
687  hEnergy ->Fill (Energy);
688  hDepth ->Fill (dep);
689  hDirCosineX ->Fill (cx);
690  hDirCosineY ->Fill (cy);
691  hDirCosineZ ->Fill (cz);
692  hTheta ->Fill (theta);
693  hPhi ->Fill (phi);
694  */
695  /*
696  // Write event by event outsputs.....
697  std::cout << "Theta: " << theta << " Phi: " << phi << " Energy: " << Energy << " Momentum: " << Momentum << std::endl;
698  std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
699  std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
700  std::cout << "Normalised..." << cx << " " << cy << " " << cz << std::endl;
701  */
702  fTree->Fill();
703  }
704 
706  // initialization
708  void MUSUN::initialization(double theta1,
709  double theta2,
710  double phi1,
711  double phi2,
712  int figflag,
713  double s_hor,
714  double s_ver1,
715  double s_ver2,
716  double& FI)
717  {
718  //
719  // Read in the data files
720  //
721  int lineNumber = 0, index = 0;
722  char inputLine[10000];
723  std::string fROOTfile;
724 
725  for (int a = 0; a < 121; ++a)
726  for (int b = 0; b < 62; ++b)
727  for (int c = 0; c < 50; ++c)
728  spmu[a][b][c] = 0;
729  for (int a = 0; a < 23401; ++a)
730  fnmu[a] = 0;
731  for (int a = 0; a < 360; ++a)
732  for (int b = 0; b < 91; ++b)
733  depth[a][b] = 0;
734  for (int a = 0; a < 360; ++a)
735  for (int b = 0; b < 91; ++b)
736  fmu[a][b] = 0;
737 
738  std::ostringstream File1LocStream;
739  File1LocStream << fInputDir << fInputFile1;
740  std::string File1Loc = File1LocStream.str();
741  cet::search_path sp1("FW_SEARCH_PATH");
742  if (sp1.find_file(fInputFile1, fROOTfile)) File1Loc = fROOTfile;
743  std::ifstream file1(File1Loc.c_str(), std::ios::in);
744  if (!file1.good())
745  throw cet::exception("MUSUNGen")
746  << "\nFile1 " << fInputFile1 << " not found in FW_SEARCH_PATH or at " << fInputDir
747  << "\n\n";
748 
749  while (file1.good()) {
750  //std::cout << "Looking at file 1...." << std::endl;
751  file1.getline(inputLine, 9999);
752  char* token;
753  token = strtok(inputLine, " ");
754  while (token != NULL) {
755  //std::cout << "While loop file 1..." << std::endl;
756  fmu[index][lineNumber] = atof(token);
757  token = strtok(NULL, " ");
758  index++;
759  if (index == 360) {
760  //std::cout << "If statement file 1..." << std::endl;
761  index = 0;
762  lineNumber++;
763  }
764  }
765  }
766  file1.close();
767 
768  std::ostringstream File2LocStream;
769  File2LocStream << fInputDir << fInputFile2;
770  std::string File2Loc = File2LocStream.str();
771  cet::search_path sp2("FW_SEARCH_PATH");
772  if (sp2.find_file(fInputFile2, fROOTfile)) File2Loc = fROOTfile;
773  std::ifstream file2(File2Loc.c_str(), std::ios::binary | std::ios::in);
774  if (!file2.good())
775  throw cet::exception("MUSUNGen")
776  << "\nFile2 " << fInputFile2 << " not found in FW_SEARCH_PATH or at " << fInputDir
777  << "\n\n";
778 
779  int i1 = 0, i2 = 0, i3 = 0;
780  float readVal;
781  while (file2.good()) {
782  //std::cout << "Looking at file 2...." << std::endl;
783  file2.read((char*)(&readVal), sizeof(float));
784  spmu[i1][i2][i3] = readVal;
785  i1++;
786  if (i1 == 121) {
787  //std::cout << "First if statement file 2..." << std::endl;
788  i2++;
789  i1 = 0;
790  }
791  if (i2 == 62) {
792  //std::cout << "Second if statement file 2..." << std::endl;
793  i3++;
794  i2 = 0;
795  }
796  }
797  file2.close();
798  for (int i = 0; i < 120; i++)
799  for (int j = 0; j < 62; j++)
800  for (int k = 0; k < 51; k++)
801  spmu[i][j][k] = spmu[i + 1][j][k];
802  spmu[1][1][0] = 0.000853544;
803  //std::cout << "Set spmu to some value..." << std::endl;
804 
805  std::ostringstream File3LocStream;
806  File3LocStream << fInputDir << fInputFile3;
807  std::string File3Loc = File3LocStream.str();
808  cet::search_path sp3("FW_SEARCH_PATH");
809  if (sp3.find_file(fInputFile3, fROOTfile)) File3Loc = fROOTfile;
810  std::ifstream file3(File3Loc.c_str(), std::ios::in);
811  if (!file3.good())
812  throw cet::exception("MUSUNGen")
813  << "\nFile3 " << fInputFile3 << " not found in FW_SEARCH_PATH or at " << fInputDir
814  << "\n\n";
815 
816  lineNumber = index = 0;
817  while (file3.good()) {
818  //std::cout << "Looking at file 3...." << std::endl;
819  file3.getline(inputLine, 9999);
820  char* token;
821  token = strtok(inputLine, " ");
822  while (token != NULL) {
823  //std::cout << "While loop file 3..." << std::endl;
824  depth[index][lineNumber] = atof(token);
825  token = strtok(NULL, " ");
826  index++;
827  if (index == 360) {
828  //std::cout << "If statement file 3..." << std::endl;
829  index = 0;
830  lineNumber++;
831  }
832  }
833  }
834  file3.close();
835 
836  //
837  // Set up variables
838  //
839 
840  the1 = theta1;
841  the2 = theta2;
842  // for c2: c1 and c2 are unused
843  //double c1 = cos(M_PI/180.*theta1);
844  //double c2 = cos(M_PI/180.*theta2);
845  ph1 = M_PI / 180. * phi1;
846  ph2 = M_PI / 180. * phi2;
847  // for c2: dph is unused
848  //double dph = ph2-ph1;
849 
850  int ipc = 1;
851  double theta = theta1;
852  double dc = 1.;
853  double sc = 0.;
854  int iteration = 0;
855  while (theta < theta2 - dc / 2.) {
856  theta += dc / 2.;
857  double theta0 = M_PI / 180. * theta;
858  double cc = cos(theta0);
859  double ash = s_hor * cc;
860  double asv01 = s_ver1 * sqrt(1. - cc * cc);
861  double asv02 = s_ver2 * sqrt(1. - cc * cc);
862  int ic1 = (theta + 0.999);
863  int ic2 = ic1 + 1;
864  if (ic2 > 91) ic2 = 91;
865  if (ic1 < 1) ic1 = 1;
866  double phi = phi1;
867  double dp = 1.;
868 
869  while (phi < phi2 - dp / 2.) {
870  phi += dp / 2.;
871  // the long side of the cavern is pointing at 14 deg to the north:
872  // double phi0 = M_PI / 180. * (phi + 90. - 14);
873 
874  // Want our co-ord system going from East to South.
875  double phi0 = M_PI / 180. * (phi + fCavernAngle);
876 
877  double asv1 = asv01 * fabs(cos(phi0));
878  double asv2 = asv02 * fabs(sin(phi0));
879  double asv0 = ash + asv1 + asv2;
880  double fl = 1.;
881  if (figflag == 1) fl = asv0;
882  int ip1 = (phi + 0.999);
883  int ip2 = ip1 + 1;
884  if (ip2 > 360) ip2 = 1;
885  if (ip1 < 1) ip1 = 360;
886  double sp1 = 0.;
887 
888  for (int ii = 0; ii < 4; ii++) {
889  int iic = ii / 2;
890  int iip = ii % 2;
891  if (ip1 == 360 && (ii == 1 || ii == 3)) iip = -359;
892  if (fmu[ip1 + iip - 1][ic1 + iic - 1] < 0) {
893  if (pow(10., fmu[ip1 + iip - 1][ic1 + iic - 1]) / 4 > 1e-6) {
894  std::cout << "Looking at fmu [ " << ip1 << " + " << iip << " - 1 (" << ip1 + iip - 1
895  << ") ] [ " << ic1 << " + " << iic << " - 1 (" << ic1 + iic - 1 << ") ] ."
896  << "\nChanging sp1 from " << sp1 << " to "
897  << sp1 + pow(10., fmu[ip1 + iip - 1][ic1 + iic - 1]) / 4 << "..........."
898  << sp1 << " + 10 ^ (" << fmu[ip1 + iip - 1][ic1 + iic - 1] << ") / 4 "
899  << std::endl;
900  }
901  sp1 = sp1 + pow(10., fmu[ip1 + iip - 1][ic1 + iic - 1]) / 4;
902  }
903  }
904  /*
905  std::cout << iteration<< " time of new sc value! Theta " << theta << ", phi " << phi + dp / 2. << ", sc = " << sc + sp1 * fl * dp * M_PI / 180. * sin(theta0) * dc * M_PI / 180. << " = "
906  << sc << " + " << sp1 << " * " << fl << " * " << dp << " * " << M_PI/180 << " * sin(" << theta0 << ") * " << dc << " * " << M_PI/180 << ".....sin(theta)=" << sin(theta) << "\n"
907  << std::endl; */
908  sc = sc + sp1 * fl * dp * M_PI / 180. * sin(theta0) * dc * M_PI / 180.;
909  ++iteration;
910  ipc = ipc + 1;
911  fnmu[ipc - 1] = sc;
912  phi = phi + dp / 2.;
913  }
914 
915  theta = theta + dc / 2.;
916  }
917  //std::cout << *FI << " = " << sc << std::endl;
918  FI = sc;
919  for (int ipc1 = 0; ipc1 < ipc; ipc1++)
920  fnmu[ipc1] = fnmu[ipc1] / fnmu[ipc - 1];
921  }
922 
924  // sampling
926  void MUSUN::sampling(double& E,
927  double& theta,
928  double& phi,
929  double& dep,
930  CLHEP::HepRandomEngine& engine)
931  {
932  CLHEP::RandFlat flat(engine);
933  CLHEP::RandGaussQ gauss(engine);
934 
935 #if 0 // this code is disabled for good
936  double xfl = flat.fire();
937  int loIndex = 0, hiIndex = 32400;
938  int i = (loIndex+hiIndex)/2;
939  bool foundIndex = false;
940  if( xfl < fnmu[loIndex] ) {
941  i = loIndex;
942  foundIndex = true;
943  } else if ( xfl > fnmu[hiIndex] ) {
944  i = hiIndex;
945  foundIndex = true;
946  } else if ( xfl > fnmu[i-1] && xfl <= fnmu[i] )
947  foundIndex = true;
948  while( !foundIndex ) {
949  if( xfl < fnmu[i] )
950  hiIndex = i;
951  else
952  loIndex = i;
953  i = (loIndex + hiIndex)/2;
954 
955  if( xfl > fnmu[i-1] && xfl <= fnmu[i] )
956  foundIndex = true;
957  }
958 #else
959  double xfl = flat.fire();
960  int i = 0;
961  while (xfl > fnmu[i])
962  ++i;
963 #endif
964  int ic = (i - 2) / 360;
965  int ip = i - 2 - ic * 360;
966 
967  xfl = flat.fire();
968  theta = the1 + ((double)ic + xfl);
969  xfl = flat.fire();
970  phi = ph1 + ((double)ip + xfl);
971  if (phi > 360) phi = phi - 360;
972  dep = depth[ip][ic] * fRockDensity;
973 
974  int ic1 = cos(M_PI / 180. * theta) * 50.;
975  if (ic1 < 0) ic1 = 0;
976  if (ic1 > 50) ic1 = 50;
977  int ip1 = dep / 200. - 16;
978  if (ip1 < 0) ip1 = 0;
979  if (ip1 > 61) ip1 = 61;
980 
981  xfl = flat.fire();
982 #if 0
983  loIndex = 0, hiIndex = 120;
984  int j = (loIndex+hiIndex)/2;
985  foundIndex = false;
986  if( xfl < spmu[loIndex][ip1][ic1] ) {
987  j = loIndex;
988  foundIndex = true;
989  } else if ( xfl > spmu[hiIndex][ip1][ic1] ) {
990  j = hiIndex;
991  foundIndex = true;
992  } else if ( xfl > spmu[j-1][ip1][ic1] && xfl <= spmu[j][ip1][ic1] )
993  foundIndex = true;
994  while( !foundIndex ) {
995  if( xfl < spmu[j][ip1][ic1] )
996  hiIndex = j;
997  else
998  loIndex = j;
999  j = (loIndex + hiIndex)/2;
1000 
1001  if( xfl > spmu[j-1][ip1][ic1] && xfl <= spmu[j][ip1][ic1] )
1002  foundIndex = true;
1003  }
1004 #else
1005  int j = 0;
1006  while (xfl > spmu[j][ip1][ic1])
1007  ++j;
1008 #endif
1009 
1010  double En1 = 0.05 * (j - 1);
1011  double En2 = 0.05 * (j);
1012  E = pow(10., En1 + (En2 - En1) * flat.fire());
1013 
1014  return;
1015  }
1016 
1017 } //end namespace evgen
1018 
1019 namespace evgen {
1020 
1022 
1023 } //end namespace evgen
base_engine_t & createEngine(seed_t seed)
unsigned int NEvents
double fEmax
Maximum Kinetic Energy (GeV)
void beginRun(art::Run &run)
double fZmax
Maximum Z position.
double fYmax
Maximum Y position.
double fDetRotZ
How much the detector should be rotated around Z axis (rotates the coordinate system clockwise)) ...
double fCavernAngle
Angle of the detector from the North to the East.
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:258
constexpr auto fullRun()
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
double fEmin
Minimum Kinetic Energy (GeV)
double fThetamin
Minimum theta.
std::string fInputDir
Input Directory.
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Run.h:121
Particle class.
static const int kGAUS
double fDetRotY
How much the detector should be rotated around Y axis (rotates the coordinate system clockwise)) ...
double fT0
Central t position (ns) in world coordinates.
Definition: Run.h:37
double fDetRotX
How much the detector should be rotated around X axis (rotates the coordinate system clockwise)) ...
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
int fTDist
How to distribute t (gaus, or uniform)
double fRockDensity
TString part[npart]
Definition: Style.C:32
int fPDG
PDG code of particles to generate.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
void endRun(art::Run &run)
void beginJob()
Definition: Breakpoints.cc:14
double fmu[360][91]
Float_t E
Definition: plot.C:20
double fPhimax
Maximum phi.
single particles thrown at the detector
Definition: MCTruth.h:26
std::string fInputFile3
Input File 3.
double fPhimin
Minimum phi.
double fChargeRatio
Charge ratio of particle / anti-particle.
std::string fInputFile2
Input File 2.
double depth[360][91]
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
ART objects.
int figflag
If want sampled from sphere or parallelepiped.
double spmu[121][62][51]
module to produce single or multiple specified particles in the detector
void initialization(double theta1, double theta2, double phi1, double phi2, int figflag, double s_hor, double s_ver1, double s_ver2, double &FI)
void SampleOne(unsigned int i, simb::MCTruth &mct, CLHEP::HepRandomEngine &engine)
double fYmin
Minimum Y position.
ifstream in
Definition: comparison.C:7
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
double fZmin
Minimum Z position.
double fSigmaT
Variation in t position (ns)
Float_t sc
Definition: plot.C:23
double fnmu[32401]
#define MF_LOG_DEBUG(id)
Definition: MVAAlg.h:12
CLHEP::HepRandomEngine & fEngine
art-managed random-number engine
double fThetamax
Maximum theta.
double fXmin
Minimum X position.
std::string const & DetectorName() const
Returns a string with the name of the detector, as configured.
Definition: GeometryCore.h:203
void produce(art::Event &evt)
TCEvent evt
Definition: DataStructs.cxx:8
Event generator information.
Definition: MCTruth.h:32
#define MF_LOG_WARNING(category)
Float_t e
Definition: plot.C:35
Namespace collecting geometry-related classes utilities.
void sampling(double &E, double &theta, double &phi, double &dep, CLHEP::HepRandomEngine &engine)
Event Generation using GENIE, cosmics or single particles.
std::string fInputFile1
Input File 1.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
double fXmax
Maximum X position.