LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
LightSource_module.cc
Go to the documentation of this file.
1 
58 #ifndef EVGEN_LIGHTSOURCE_H
59 #define EVGEN_LIGHTSOURCE_H
60 
61 // C++ includes.
62 #include <iostream>
63 #include <string>
64 #include <cmath>
65 #include <memory>
66 #include <fstream>
67 
68 // ART includes
71 #include "fhiclcpp/ParameterSet.h"
79 #include "cetlib_except/exception.h"
80 
81 // art extensions
83 
84 // nutools includes
88 
89 // lar includes
95 
96 #include "TVector3.h"
97 #include "TLorentzVector.h"
98 #include "TTree.h"
99 
100 #include "CLHEP/Random/RandFlat.h"
101 #include "CLHEP/Random/RandGaussQ.h"
102 
103 namespace evgen {
104 
106  class LightSource : public art::EDProducer {
107  public:
108  explicit LightSource(fhicl::ParameterSet const& pset);
109 
110  void produce(art::Event & evt);
111  void beginRun(art::Run& run);
112 
113  private:
114 
115  void Sample(simb::MCTruth &truth);
116 
117  // for c2: fSeed is unused
118  //int fSeed; //random number seed
119  std::string fVersion; //version of the configuration
120 
121  // Flags to mark module modes
122  static const int kUNIF = 0;
123  static const int kGAUS = 1;
124  static const int kFILE = 0;
125  static const int kSCAN = 1;
126 
127  // File stream, filename and empty string for file processing
128  std::ifstream fInputFile;
129  std::string fFileName;
130  char fDummyString[256];
131 
132  // A ttree to keep track of where particles have been shot - ends up in histos.root
134  TLorentzVector fShotPos;
135  TLorentzVector fShotMom;
136  Int_t fEvID;
137 
138  // Parameters loaded from config - both modes
139  int fSourceMode; // Mode to run in - scan or file
140  bool fFillTree; // Do we want to create a TTree of shot particles?
141  int fPosDist; //
142  int fTDist; // Random distributions to use : 1= gauss, 0= uniform
143  int fPDist; //
144 
145  //Scan mode specific parameters
146  int fXSteps; //
147  int fYSteps; // Number of steps to take in each dimension
148  int fZSteps; //
149 
150  sim::PhotonVoxelDef fThePhotonVoxelDef; // The photon voxel definition object for scan mode
151 
152  int fVoxelCount; // Total number of voxels
153  int fCurrentVoxel; // Counter to keep track of vox ID
154 
155 
156  // TPC Measurements
157  TVector3 fTPCCenter;
158  TVector3 fTPCDimensions;
159  std::vector<double> fRegionMin;
160  std::vector<double> fRegionMax;
162 
163 
164  // Parameters used to shoot in distributions
165  double fX; // central x position of source
166  double fY; // central y position of source
167  double fZ; // central z position of source
168  bool fPointSource; // Point-like light source in fX, fY, fZ
169  double fT; // central t position of source
170  double fSigmaX; // x width
171  double fSigmaY; // y width
172  double fSigmaZ; // z width
173  double fSigmaT; // t width
174  double fP; // central momentm of photon
175  double fSigmaP; // mom width;
176 
177  // Number of photons per event
178  int fN; // number of photons per event
179 
182  };
183 }
184 
185 namespace evgen{
186 
187  //----------------------------------------------------------------
189  {
190 
191  // get the random number seed, use a random default if not specified
192  // in the configuration file.
193 
194  fSourceMode = (pset.get<int >("SourceMode") );
195  fFillTree = (pset.get<bool>("FillTree") );
196  fPosDist = (pset.get<int >("PosDist") );
197  fPDist = (pset.get<int >("PDist") );
198  fTDist = (pset.get<int >("TDist") );
199 
200  // create a default random engine; obtain the random seed from NuRandomService,
201  // unless overridden in configuration with key "Seed"
203  ->createEngine(*this, pset, "Seed");
204 
205  // load optional parameters in function
206  produces< sumdata::RunData, art::InRun >();
207  produces< std::vector<simb::MCTruth> >();
208 
209  if(fSourceMode==kFILE)
210  {
211  fFileName = pset.get<std::string>("SteeringFile");
212  fInputFile.open(fFileName.c_str());
213  fInputFile.getline(fDummyString,256);
214 
215  }
216  else if (fSourceMode==kSCAN)
217  {
218  fT = pset.get<double>("T0");
219  fSigmaT = pset.get<double>("SigmaT");
220  fN = pset.get<int >("N");
221 
222  fFirstVoxel = pset.get<int >("FirstVoxel");
223  fLastVoxel = pset.get<int >("LastVoxel");
224 
225  fP = pset.get<double>("P");
226  fSigmaP = pset.get<double>("SigmaP");
227 
228  fUseCustomRegion = pset.get<bool>("UseCustomRegion");
229  fPointSource = pset.get<bool>("PointSource",false);
230 
231  if(fUseCustomRegion)
232  {
233  fRegionMin = pset.get< std::vector<double> >("RegionMin");
234  fRegionMax = pset.get< std::vector<double> >("RegionMax");
235  fXSteps = pset.get<int >("XSteps");
236  fYSteps = pset.get<int >("YSteps");
237  fZSteps = pset.get<int >("ZSteps");
238  }
239 
241  // get TPC dimensions removed. -TA
242 
243 
244  fCurrentVoxel=0;
245 
246  // define voxelization based on parameters read from config.
247  // There are two modes - either read the dimensions of the TPC from
248  // the geometry, or use values specified by the user.
249  if(!fUseCustomRegion)
250  {
253  }
254  else
255  {
257  fRegionMax[0],
258  fXSteps,
259  fRegionMin[1],
260  fRegionMax[1],
261  fYSteps,
262  fRegionMin[2],
263  fRegionMax[2],
264  fZSteps);
265  }
266 
267 
268  // Set distribution widths to voxel size
269 
273 
274  // Get number of voxels we will step through
275 
277 
278  if(fLastVoxel<0) fLastVoxel = fVoxelCount;
279 
280  mf::LogVerbatim("LightSource") << "Light Source : Determining voxel params : "
281  << fVoxelCount << " "
282  << fSigmaX << " "
283  << fSigmaY << " "
284  <<fSigmaZ;
285 
286  }
287  else{
288  throw cet::exception("LightSource") << "EVGEN Light Source : Unrecognised light source mode\n";
289  }
290 
291 
292 
293 
294  if(fFillTree)
295  {
297  fPhotonsGenerated = tfs->make<TTree>("PhotonsGenerated","PhotonsGenerated");
298  fPhotonsGenerated->Branch("X",&(fShotPos[0]),"X/D");
299  fPhotonsGenerated->Branch("Y",&(fShotPos[1]),"Y/D");
300  fPhotonsGenerated->Branch("Z",&(fShotPos[2]),"Z/D");
301  fPhotonsGenerated->Branch("T",&(fShotPos[3]),"T/D");
302  fPhotonsGenerated->Branch("PX",&(fShotMom[0]),"PX/D");
303  fPhotonsGenerated->Branch("PY",&(fShotMom[1]),"PY/D");
304  fPhotonsGenerated->Branch("PZ",&(fShotMom[2]),"PZ/D");
305  fPhotonsGenerated->Branch("PT",&(fShotMom[3]),"PT/D");
306  fPhotonsGenerated->Branch("EventID",&fEvID,"EventID/I");
307  }
308  }
309 
310 
311  //____________________________________________________________________________
313  {
314 
315  // grab the geometry object to see what geometry we are using
317  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
318 
319  run.put(std::move(runcol));
320 
322 
323  return;
324  }
325 
326  //----------------------------------------------------------------
328  {
329 
330  // FILE MODE -
331  // Each event, read coordinates of gun and number of photons to shoot from file
332 
333  if(fSourceMode==kFILE)
334  {
335  // Loop file if required
336  if(fInputFile.eof()){
337  mf::LogWarning("LightSource") << "EVGEN Light Source : Warning, reached end of file,"
338  << " looping back to beginning";
339  fInputFile.seekg(0,std::ios::beg);
340  fInputFile.clear();
341  }
342 
343  if(!fInputFile.is_open() || fInputFile.fail() ){
344  throw cet::exception("LightSource") << "EVGEN Light Source : File error in "
345  << fFileName << "\n";
346  }
347  else{
348  // read in one line
349  fInputFile >> fX >> fY >> fZ >> fT
350  >> fSigmaX >> fSigmaY >> fSigmaZ >> fSigmaT
351  >> fP >> fSigmaP >> fN;
352  fInputFile.getline(fDummyString,256);
354  fX + fSigmaX,
355  1,
356  fY - fSigmaY,
357  fY + fSigmaY,
358  1,
359  fZ - fSigmaZ,
360  fZ + fSigmaZ,
361  1);
362 
363  fCurrentVoxel=0;
364  }
365  }
366 
367 
368  // SCAN MODE -
369  // Step through detector using a number of steps provided in the config file
370  // firing a constant number of photons from each point
371  else if(fSourceMode==kSCAN)
372  {
373 
374 
375  TVector3 VoxelCenter = fThePhotonVoxelDef.GetPhotonVoxel(fCurrentVoxel).GetCenter();
376 
377  fX = VoxelCenter.X();
378  fY = VoxelCenter.Y();
379  fZ = VoxelCenter.Z();
380 
381  }
382 
383 
384  // UNRECOGNISED MODE
385  // - neither file or scan mode, probably a config file error
386 
387  else{
388  throw cet::exception("LightSource") <<"EVGEN : Light Source, unrecognised source mode\n";
389  }
390 
391 
392 
393  // std::cout<<"EVGEN Light source to be placed at (x, y, z, t, dx, dy, dz, dt, n) " <<
394  // fX << " " << fY << " " <<fZ << " " << fT << " " <<
395  // fSigmaX << " " << fSigmaY << " " << fSigmaZ << " " << fSigmaT << " " <<
396  // fN << std::endl<<std::endl;
397 
398 
399  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
400 
401  simb::MCTruth truth;
403  Sample(truth);
404 
405  // std::cout << "put mctruth into the vector" << std::endl;
406  truthcol->push_back(truth);
407 
408  // std::cout << "add vector to the event " << truthcol->size() << std::endl;
409  evt.put(std::move(truthcol));
410 
411  phot::PhotonVisibilityService* vis = nullptr;
412  try {
414  }
415  catch (art::Exception const& e) {
416  // if the service is not configured, then this is not a build job
417  // (it is a simple generation job instead)
418  if (e.categoryCode() != art::errors::ServiceNotFound) throw;
419  }
420 
421  if(vis && vis->IsBuildJob())
422  {
423  mf::LogVerbatim("LightSource") << "Light source : Stowing voxel params ";
425  }
426 
428  {
429  ++fCurrentVoxel;
430  }
431  else
432  {
433  mf::LogVerbatim("LightSource") << "EVGEN Light Source fully scanned detector. Starting over.";
435  }
436 
437  return;
438 
439 
440  }
441 
442 
444  {
445  mf::LogVerbatim("LightSource") <<"Light source debug : Shooting at " << fX <<" " << fY<<" "<< fZ;
446 
447  // get the random number generator service and make some CLHEP generators
449  CLHEP::HepRandomEngine &engine = rng->getEngine();
450  CLHEP::RandFlat flat(engine);
451  CLHEP::RandGaussQ gauss(engine);
452 
453  for(int j=0; j!=fN; ++j){
454  // Choose momentum (supplied in eV, convert to GeV)
455  double p = fP;
456  if (fPDist == kGAUS) {
457  p = gauss.fire(fP, fSigmaP);
458  }
459  else {
460  p = fP + fSigmaP*(2.0*flat.fire()-1.0);
461  }
462  p /= 1000000000.;
463 
464  // Choose position
465  TVector3 x;
466  if(fPointSource) {
467  x[0] = fX;
468  x[1] = fY;
469  x[2] = fZ;
470  }
471  else {
472  if (fPosDist == kGAUS) {
473  x[0] = gauss.fire(fX, fSigmaX);
474  x[1] = gauss.fire(fY, fSigmaY);
475  x[2] = gauss.fire(fZ, fSigmaZ);
476  }
477  else {
478  x[0] = fX + fSigmaX*(2.0*flat.fire()-1.0);
479  x[1] = fY + fSigmaY*(2.0*flat.fire()-1.0);
480  x[2] = fZ + fSigmaZ*(2.0*flat.fire()-1.0);
481  }
482 
483  }
484 
485  // Choose time
486  double t;
487  if (fTDist == kGAUS) {
488  t = gauss.fire(fT, fSigmaT);
489  }
490  else {
491  t = fT + fSigmaT * (2.0 * flat.fire()-1.0);
492  }
493 
494 
495  //assume the position is relative to the center of the TPC
496  //x += fTPCCenter;
497 
498  fShotPos = TLorentzVector(x[0], x[1], x[2], t);
499 
500 
501  // Choose angles
502  double costh = 2 * flat.fire() - 1;
503  double sinth = pow(1-pow(costh,2),0.5);
504  double phi = 2 * M_PI * flat.fire();
505 
506  // Generate momentum 4-vector
507 
508  fShotMom = TLorentzVector( p*sinth*cos(phi),
509  p*sinth*sin(phi),
510  p*costh,
511  p );
512 
513  int trackid = -1*(j+1); // set track id to -i as these are all primary particles and have id <= 0
514  std::string primary("primary");
515  int PDG=0; //optical photons have PDG 0
516 
517  simb::MCParticle part(trackid, PDG, primary);
519 
520  if(fFillTree)
521  fPhotonsGenerated->Fill();
522 
523  mct.Add(part);
524  }
525 
526  }
527 
528 
529 }
530 
531 
532 namespace evgen{
533 
535 
536 }//end namespace evgen
537 
538 
539 #endif // EVGEN_LIGHTSOURCE_H
540 
542 
Float_t x
Definition: compare.C:6
std::vector< double > fRegionMin
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
LightSource(fhicl::ParameterSet const &pset)
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:261
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:78
static const int kSCAN
static const int kGAUS
std::vector< double > fRegionMax
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
TVector3 GetCenter() const
Particle class.
void Add(simb::MCParticle &part)
Definition: MCTruth.h:77
Definition: Run.h:30
void produce(art::Event &evt)
int GetNVoxels() const
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
base_engine_t & getEngine() const
base_engine_t & createEngine(seed_t seed)
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
single particles thrown at the detector
Definition: MCTruth.h:24
T get(std::string const &key) const
Definition: ParameterSet.h:231
TString part[npart]
Definition: Style.C:32
TVector3 GetVoxelSize() const
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
void beginRun(art::Run &run)
sim::PhotonVoxelDef fThePhotonVoxelDef
void Sample(simb::MCTruth &truth)
static const int kUNIF
const sim::PhotonVoxelDef & GetVoxelDef() const
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
T * make(ARGS...args) const
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
static const int kFILE
A module for optical MC testing and library building.
TCEvent evt
Definition: DataStructs.cxx:5
Event generator information.
Definition: MCTruth.h:30
Float_t e
Definition: plot.C:34
Namespace collecting geometry-related classes utilities.
Event Generation using GENIE, cosmics or single particles.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
PhotonVoxel GetPhotonVoxel(int ID) const