LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
FileMuons_module.cc
Go to the documentation of this file.
1 // C++ includes.
10 #include <string>
11 #include <cmath>
12 #include <memory>
13 
14 // Framework includes
16 #include "fhiclcpp/ParameterSet.h"
23 
24 // nutools includes
27 
28 // lar includes
31 
32 #include "TVector3.h"
33 #include "TDatabasePDG.h"
34 
36 
37 
38 #include <vector>
39 #include <string>
41 #include "TFile.h"
42 #include "TTree.h"
43 
44 #include <stdio.h>
45 #include <iostream>
46 #include <fstream>
47 
48 namespace simb { class MCTruth; }
49 
50 namespace evgen {
51 
53  class FileMuons : public art::EDProducer {
54 
55  public:
56  explicit FileMuons(fhicl::ParameterSet const& pset);
57  virtual ~FileMuons();
58 
59  // This is called for each event.
60  void produce(art::Event& evt);
61  void beginJob();
62  void beginRun(art::Run& run);
63  void endJob();
64 
65  private:
66 
67  void ReadEvents(simb::MCTruth &mct);
68 
69  int fEventNumberOffset; // Where in file to start.
70  std::vector<int> fPDG;
71  std::vector<double> fXYZ_Off;
72  std::string fFileName;
73  std::string fMuonsFileType;
74  std::string fTreeName;
75  std::vector<std::string> fBranchNames;
76 
77  std::ifstream *fMuonFile;
78  TFile *fMuonFileR;
79  TTree *TNtuple;
80  unsigned int countFile;
81 
82  Float_t xtmp, ytmp, ztmp;
83  Float_t pxtmp, pytmp, pztmp;
84  Float_t charge;
85  Float_t E;
86  Float_t costheta;
87  Float_t phi;
88  Float_t xdet;
89  Float_t ydet;
90  Float_t zdet;
91 
92  TBranch *b_x;
93  TBranch *b_y;
94  TBranch *b_z;
95  TBranch *b_E;
96  TBranch *b_costheta;
97  TBranch *b_phi;
98  TBranch *b_xdet;
99  TBranch *b_ydet;
100  TBranch *b_zdet;
101  TBranch *b_px;
102  TBranch *b_py;
103  TBranch *b_pz;
104  TBranch *b_charge;
105 
106 
107  };
108 } // namespace
109 
110 
111 
112 namespace evgen{
113 
114  //____________________________________________________________________________
115  FileMuons::FileMuons(fhicl::ParameterSet const& pset)
116  : fEventNumberOffset(pset.get< int >("EventNumberOffset"))
117  , fPDG (pset.get< std::vector<int> >("PDG") )
118  , fXYZ_Off (pset.get< std::vector<double> >("InitialXYZOffsets"))
119  , fFileName (pset.get< std::string >("FileName") )
120  , fMuonsFileType (pset.get< std::string >("MuonsFileType") )
121  , fTreeName (pset.get< std::string >("TreeName") )
122  , fBranchNames (pset.get< std::vector<std::string> >("BranchNames") )
123  {
124 
125  produces< std::vector<simb::MCTruth> >();
126  produces< sumdata::RunData, art::InRun >();
127 
128  }
129 
130  //____________________________________________________________________________
132  {
134  mf::LogInfo("FileMuons : starting at event ") << countFile <<std::endl;
135 
136  if (fMuonsFileType.compare("source")==0)
137  {
138  std::cout << "FileMuons: Not yet equipped to walk through muons with TFS mojo."<< std::endl;
139  }
140  else if (fMuonsFileType.compare("root")==0)
141  {
142  std::cout << "FileMuons: You have chosen to read muons from Root File " << fFileName << std::endl;
143  }
144  else if (fMuonsFileType.compare("text")==0)
145  {
146  std::cout << "FileMuons: You have chosen to read muons from " << fFileName << "." << std::endl;
147  }
148  else
149  {
150  std::cout << "FileMuons: You must specify one of source/text/root file to read for muons."<< std::endl;
151 
152  }
153 
154  if (fMuonsFileType.compare("text")==0)
155  {
156  fMuonFile = new std::ifstream(fFileName.c_str());
157  long begin = fMuonFile->tellg();
158  fMuonFile->seekg (0, std::ios::end);
159  long end = fMuonFile->tellg();
160  std::cout << "FileMuons: " << fFileName << " size is: " << (end-begin) << " bytes.\n";
161  fMuonFile->seekg (0, std::ios::beg);
162 
163  for (unsigned int header=0; header<3 && fMuonFile->good(); ++header)
164  {
165  std::string line;
166  getline(*fMuonFile,line);
167  }
168  if (!fMuonFile->good())
169  {
170  std::cout << "FileMuons: Problem reading muon file header."<< std::endl;
171  }
172  } // fMuonsFileType is a text file.
173  else if (fMuonsFileType.compare("root")==0)
174  {
175  fMuonFileR = new TFile(fFileName.c_str(),"READ");
176  TNtuple = (TTree*)(fMuonFileR->Get(fTreeName.c_str()));
177 
178  TNtuple->SetBranchAddress("x", &xtmp, &b_x);
179  TNtuple->SetBranchAddress("y", &ytmp, &b_y);
180  TNtuple->SetBranchAddress("z", &ztmp, &b_z);
181  TNtuple->SetBranchAddress("E", &E, &b_E);
182  TNtuple->SetBranchAddress("costheta", &costheta, &b_costheta);
183  TNtuple->SetBranchAddress("phi", &phi, &b_phi);
184  TNtuple->SetBranchAddress("xdet", &xdet, &b_xdet);
185  TNtuple->SetBranchAddress("ydet", &ydet, &b_ydet);
186  TNtuple->SetBranchAddress("zdet", &zdet, &b_zdet);
187  TNtuple->SetBranchAddress("px", &pxtmp, &b_px);
188  TNtuple->SetBranchAddress("py", &pytmp, &b_py);
189  TNtuple->SetBranchAddress("pz", &pztmp, &b_pz);
190  TNtuple->SetBranchAddress("charge", &charge, &b_charge);
191 
192  } // fMuonsFileType is a root file.
193 
194  }
195 
196  //____________________________________________________________________________
198  {
199  if (fMuonsFileType.compare("text")==0)
200  fMuonFile->close();
201  if (fMuonsFileType.compare("root")==0)
202  fMuonFileR->Close();
203  }
204 
205  //____________________________________________________________________________
207  {
208 
209  // grab the geometry object to see what geometry we are using
211  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
212 
213  run.put(std::move(runcol));
214 
215  return;
216  }
217 
218  //____________________________________________________________________________
220  {
221  }
222 
223  //____________________________________________________________________________
225  {
226 
228  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
229 
230  simb::MCTruth truth;
232  ReadEvents(truth);
233 
234 // std::cout << "put mctruth into the vector" << std::endl;
235  truthcol->push_back(truth);
236 
237 // std::cout << "add vector to the event " << truthcol->size() << std::endl;
238  evt.put(std::move(truthcol));
239 
240  return;
241  }
242 
243  //____________________________________________________________________________
245  {
246 
247 // std::cout << "size of particle vector is " << fPDG.size() << std::endl;
248 
250  for (unsigned int i=0; i<fPDG.size(); ++i) {
251 
252  // Choose momentum
253  //double p = 0.0;
254  double m(0.108);
255 
256  TVector3 x;
257  TVector3 p;
258  Double_t q = 0.;
259  Int_t pdgLocal;
260 
261  if (fMuonsFileType.compare("text")==0)
262  {
263 
264  std::string line;
265  getline(*fMuonFile,line);
266  if (!fMuonFile->good())
267  {
268  std::cout << "FileMuons: Problem reading muon file line ...."<< countFile << ". Perhaps you've exhausted the events in " << fFileName << std::endl; exit(0);
269  }
270  else
271  {
272  // std::cout << "FileMuons: getline() gives "<< line << " for event " << countFile << std::endl;
273  }
274  countFile++;
275 
276  LOG_DEBUG("FileMuons: countFile is ") << countFile <<std::endl;
277  char * cstr, *ptok;
278 
279  // Split this line into tokens
280  cstr = new char [line.size()+1];
281  strcpy (cstr, line.c_str());
282  // cstr now contains a c-string copy of str
283  ptok=strtok (cstr,"*");
284  unsigned int fieldCount = 0;
285  unsigned int posIndex = 0;
286  unsigned int pIndex = 0;
287  while (ptok!=NULL)
288  {
289  // std::cout << ptok << std::endl;
290  ptok=strtok(NULL,"*");
291  if (fieldCount==9 || fieldCount==10 || fieldCount==11)
292  {
293  p[pIndex] = atof(ptok); pIndex++;
294  // std::cout << ptok << std::endl;
295  }
296  if (fieldCount==6 || fieldCount==7 || fieldCount==8)
297  {
298  x[posIndex] = atof(ptok);
299  // make the z axis point up for x, as with p
300  if (posIndex==2) {x[posIndex] = -1.0*x[posIndex];}
301  posIndex++;
302  }
303  if (fieldCount==12)
304  {
305  q = atof(ptok);
306  }
307  fieldCount++;
308  }
309 
310  delete[] cstr;
311 
312  }
313  else if (fMuonsFileType.compare("root")==0) // from root file
314  {
315  /*
316  // Don't use this yet. Keep the specific branch-by-branch identification.
317  for (unsigned int ii=0;ii<fBranchNames.size();ii++)
318  {
319  TNtuple->SetBranchAddress(fBranchNames[ii], x+ii);
320  }
321  */
322  // TNtuple->ResetBranchAddresses();
323  TNtuple->GetEntry(countFile);
324  //TNtuple->Show(countFile);
325 
326  x.SetXYZ(xdet,ydet,-zdet); // as with txt file, make z point up.
327  // Watch for units change to mm in Modern JdJ files!!
328  // This is for pre Spring-2012 JdJ Ntuples.
329  p.SetXYZ(pxtmp,pytmp,pztmp);
330  q = charge;
331 
332  countFile++;
333 
334  } // End read.
335 
336  static TDatabasePDG pdgt;
337  pdgLocal = -q*fPDG[i];
338 
339  TParticlePDG* pdgp = pdgt.GetParticle(pdgLocal);
340  if (pdgp) m = pdgp->Mass();
341 
342 
343  // std::cout << "set the position "<<std::endl;
344  // This gives coordinates at the center of the 300mx300m plate that is 3m above top of
345  // cavern. Got these by histogramming deJong's xdet,ydet,zdet.
346  const double cryoGap = 15.0;
347  x[0] -= fXYZ_Off[0];
348  x[1] -= fXYZ_Off[1];
349  x[2] -= fXYZ_Off[2]; // 3 for plate height above top of cryostat.
350  // Now, must rotate to TPC coordinates. Let's orient TPC axis along z axis,
351  // Cosmics, mostly going along deJong's +z axis must be going along TPC -y axis.
352  x.RotateX(-M_PI/2);
353  p.RotateX(-M_PI/2);
354  //add vector of the position of the center of the point between Cryostats
355  // level with top. (To which I've added 3m - in above code - in height.)
356  // This is referenced from origin at center-right of first cryostat.
358  TVector3 off3(geom->CryostatHalfWidth()*0.01,geom->CryostatHalfHeight()*0.01,geom->CryostatLength()*0.01+cryoGap*0.01/2.0) ;
359  x += off3;
360 
361  TLorentzVector pos(x[0]*100.0, x[1]*100.0, x[2]*100.0, 0.0);
362  TLorentzVector pvec(p[0]*1000.0,p[1]*1000.0,p[2]*1000.0,std::sqrt(p.Mag2()*1000.0*1000.0+m*m));
363  std::cout << "x[m] and p [TeV] are " << std::endl;
364  x.Print();
365  p.Print();
366 
367  int trackid = -1*(i+1); // set track id to -i as these are all primary particles and have id <= 0
368  std::string primary("primary");
369  simb::MCParticle part(trackid, pdgLocal, primary);
370  part.AddTrajectoryPoint(pos, pvec);
371 
372  // std::cout << "add the particle to the primary" << std::endl;
373 
374  mct.Add(part);
375 
376  }//end loop over particles
377 
378  return;
379  }
380 
382 
383 }//end namespace evgen
Float_t x
Definition: compare.C:6
std::vector< double > fXYZ_Off
geo::Length_t CryostatHalfHeight(geo::CryostatID const &cid) const
Returns the height of the cryostat (y direction)
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:261
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:78
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::string fMuonsFileType
void ReadEvents(simb::MCTruth &mct)
void beginRun(art::Run &run)
STL namespace.
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
std::ifstream * fMuonFile
Particle class.
void Add(simb::MCParticle &part)
Definition: MCTruth.h:77
Definition: Run.h:30
std::vector< std::string > fBranchNames
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
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
unsigned int countFile
std::vector< evd::details::RawDigitInfo_t >::const_iterator begin(RawDigitCacheDataClass const &cache)
TString part[npart]
Definition: Style.C:32
std::string fFileName
module to produce single or multiple specified particles in the detector
Framework includes.
geo::Length_t CryostatHalfWidth(geo::CryostatID const &cid) const
Returns the half width of the cryostat (x direction)
#define LOG_DEBUG(id)
geo::Length_t CryostatLength(geo::CryostatID const &cid) const
Returns the length of the cryostat (z direction)
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
Event generator information.
Definition: MCTruth.h:30
Namespace collecting geometry-related classes utilities.
Event Generation using GENIE, cosmics or single particles.
std::vector< int > fPDG
art framework interface to geometry description
std::string fTreeName
void produce(art::Event &evt)