Public Member Functions | Private Member Functions | Private Attributes | List of all members
evgen::CosmicsGen Class Reference

A module to check the results from the Monte Carlo generator. More...

Inheritance diagram for evgen::CosmicsGen:
art::EDProducer art::detail::Producer art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper

Public Member Functions

 CosmicsGen (fhicl::ParameterSet const &pset)
 
- Public Member Functions inherited from art::EDProducer
 EDProducer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDProducer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Producer
virtual ~Producer () noexcept
 
 Producer (fhicl::ParameterSet const &)
 
 Producer (Producer const &)=delete
 
 Producer (Producer &&)=delete
 
Produceroperator= (Producer const &)=delete
 
Produceroperator= (Producer &&)=delete
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Modifier
 ~Modifier () noexcept
 
 Modifier ()
 
 Modifier (Modifier const &)=delete
 
 Modifier (Modifier &&)=delete
 
Modifieroperator= (Modifier const &)=delete
 
Modifieroperator= (Modifier &&)=delete
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Member Functions

void produce (art::Event &evt) override
 
void beginJob () override
 
void beginRun (art::Run &run) override
 

Private Attributes

std::vector< double > fbuffbox
 
TH2F * fPhotonAngles
 Photon rate vs angle. More...
 
TH2F * fPhotonAnglesLo
 Photon rate vs angle, low momenta. More...
 
TH2F * fPhotonAnglesMi
 Photon rate vs angle, middle momenta. More...
 
TH2F * fPhotonAnglesHi
 Photon rate vs angle, high momenta. More...
 
TH1F * fPhotonCosQ
 Photon rate vs cos(Q) More...
 
TH1F * fPhotonEnergy
 Photon energy (GeV) More...
 
TH1F * fPhotonsPerSample
 number of photons in the sampled time window More...
 
TH1F * fPhotonsInCStat
 
TH1F * fPhotonsInTPC
 
TH2F * fElectronAngles
 Electron rate vs angle. More...
 
TH2F * fElectronAnglesLo
 Electron rate vs angle, low momenta. More...
 
TH2F * fElectronAnglesMi
 Electron rate vs angle, middle momenta. More...
 
TH2F * fElectronAnglesHi
 Electron rate vs angle, high momenta. More...
 
TH1F * fElectronCosQ
 Electron rate vs cos(Q) More...
 
TH1F * fElectronEnergy
 Electron energy (GeV) More...
 
TH1F * fElectronsPerSample
 number of electrons in the sampled time window More...
 
TH1F * fElectronsInCStat
 
TH1F * fElectronsInTPC
 
TH2F * fMuonAngles
 Muon rate vs angle. More...
 
TH2F * fMuonAnglesLo
 Muon rate vs angle, low momenta. More...
 
TH2F * fMuonAnglesMi
 Muon rate vs angle, middle momenta. More...
 
TH2F * fMuonAnglesHi
 Muon rate vs angle, high momenta. More...
 
TH1F * fMuonCosQ
 Muon rate vs cos(Q) More...
 
TH1F * fMuonEnergy
 Muon energy (GeV) More...
 
TH1F * fMuonsPerSample
 number of muons in the sampled time window More...
 
TH1F * fMuonsInCStat
 
TH1F * fMuonsInTPC
 
CLHEP::HepRandomEngine & fEngine
 art-managed random-number engine More...
 
evgb::CRYHelper fCRYHelp
 CRY generator object. More...
 

Additional Inherited Members

- Public Types inherited from art::EDProducer
using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
- Public Types inherited from art::detail::Producer
template<typename UserConfig , typename KeysToIgnore = void>
using Table = Modifier::Table< UserConfig, KeysToIgnore >
 
- Public Types inherited from art::Modifier
template<typename UserConfig , typename UserKeysToIgnore = void>
using Table = ProducerTable< UserConfig, detail::ModuleConfig, UserKeysToIgnore >
 
- Static Public Member Functions inherited from art::EDProducer
static void commitEvent (EventPrincipal &ep, Event &e)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

A module to check the results from the Monte Carlo generator.

Definition at line 36 of file CosmicsGen_module.cc.

Constructor & Destructor Documentation

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

Definition at line 91 of file CosmicsGen_module.cc.

92  : art::EDProducer{pset}
93  , fbuffbox{pset.get<std::vector<double>>("BufferBox",{0.0, 0.0, 0.0, 0.0, 0.0, 0.0})}
94  // create a default random engine; obtain the random seed from NuRandomService,
95  // unless overridden in configuration with key "Seed"
96  , fEngine(art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this, pset, "Seed"))
97  , fCRYHelp{pset,
98  fEngine,
100  {
101  produces< std::vector<simb::MCTruth> >();
102  produces< sumdata::RunData, art::InRun >();
103  }
const std::string GetWorldVolumeName() const
Return the name of the world volume (needed by Geant4 simulation)
CLHEP::HepRandomEngine & fEngine
art-managed random-number engine
evgb::CRYHelper fCRYHelp
CRY generator object.
std::vector< double > fbuffbox

Member Function Documentation

void evgen::CosmicsGen::beginJob ( )
overrideprivatevirtual

Reimplemented from art::EDProducer.

Definition at line 106 of file CosmicsGen_module.cc.

107  {
109 
110  fPhotonAngles = tfs->make<TH2F>("fPhotonAngles", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
111  fPhotonAnglesLo = tfs->make<TH2F>("fPhotonAnglesLo", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
112  fPhotonAnglesMi = tfs->make<TH2F>("fPhotonAnglesMi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
113  fPhotonAnglesHi = tfs->make<TH2F>("fPhotonAnglesHi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
114  fPhotonCosQ = tfs->make<TH1F>("fPhotonCosQ", ";cos#theta;tracks", 50,-1.0,1.0);
115  fPhotonEnergy = tfs->make<TH1F>("fPhotonEnergy", ";E (GeV)", 5000,0.0,1000.0);
116  fPhotonsPerSample = tfs->make<TH1F>("fPhotonsPerSample", ";Number Photons;Samples", 100, 0, 1000);
117  fPhotonsInCStat = tfs->make<TH1F>("fPhotonsInCryostat", ";Number Photons;Samples", 100, 0, 1000);
118  fPhotonsInTPC = tfs->make<TH1F>("fPhotonsInTPC", ";Number Photons;Samples", 100, 0, 1000);
119 
120  fElectronAngles = tfs->make<TH2F>("fElectronAngles", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
121  fElectronAnglesLo = tfs->make<TH2F>("fElectronAnglesLo", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
122  fElectronAnglesMi = tfs->make<TH2F>("fElectronAnglesMi", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
123  fElectronAnglesHi = tfs->make<TH2F>("fElectronAnglesHi", ";#phi;cos#theta",36,-180.0,180.0,50,-1.0,1.0);
124  fElectronCosQ = tfs->make<TH1F>("fElectronCosQ", ";cos#theta;tracks",50,-1.0,1.0);
125  fElectronEnergy = tfs->make<TH1F>("fElectronEnergy", ";E (GeV)", 5000,0.0,1000.0);
126  fElectronsPerSample = tfs->make<TH1F>("fElectronsPerSample", ";Number Electrons;Samples", 100, 0, 1000);
127  fElectronsInCStat = tfs->make<TH1F>("fElectronsInCryotat", ";Number Electrons;Samples", 100, 0, 1000);
128  fElectronsInTPC = tfs->make<TH1F>("fElectronsInTPC", ";Number Electrons;Samples", 100, 0, 1000);
129 
130  fMuonAngles = tfs->make<TH2F>("fMuonAngles", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
131  fMuonAnglesLo = tfs->make<TH2F>("fMuonAnglesLo", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
132  fMuonAnglesMi = tfs->make<TH2F>("fMuonAnglesMi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
133  fMuonAnglesHi = tfs->make<TH2F>("fMuonAnglesHi", ";#phi;cos#theta", 36,-180.0,180.0,50,-1.0,1.0);
134  fMuonCosQ = tfs->make<TH1F>("fMuonCosQ", ";cos#theta;tracks",50,-1.0,1.0);
135  fMuonEnergy = tfs->make<TH1F>("fMuonEnergy", ";E (GeV)", 5000,0.0,1000.0);
136  fMuonsPerSample = tfs->make<TH1F>("fMuonsPerSample", ";Number Muons;Samples", 100, 0, 1000);
137  fMuonsInCStat = tfs->make<TH1F>("fMuonsInCryostat", ";Number Muons;Samples", 100, 0, 1000);
138  fMuonsInTPC = tfs->make<TH1F>("fMuonsInTPC", ";Number Muons;Samples", 100, 0, 1000);
139  }
TH2F * fElectronAnglesHi
Electron rate vs angle, high momenta.
TH2F * fPhotonAngles
Photon rate vs angle.
TH2F * fMuonAnglesHi
Muon rate vs angle, high momenta.
TH1F * fElectronEnergy
Electron energy (GeV)
TH2F * fElectronAnglesLo
Electron rate vs angle, low momenta.
TH2F * fElectronAngles
Electron rate vs angle.
TH1F * fElectronCosQ
Electron rate vs cos(Q)
TH2F * fMuonAngles
Muon rate vs angle.
TH2F * fMuonAnglesLo
Muon rate vs angle, low momenta.
TH1F * fPhotonEnergy
Photon energy (GeV)
TH2F * fElectronAnglesMi
Electron rate vs angle, middle momenta.
TH1F * fMuonCosQ
Muon rate vs cos(Q)
TH1F * fElectronsPerSample
number of electrons in the sampled time window
TH1F * fPhotonsPerSample
number of photons in the sampled time window
TH1F * fPhotonCosQ
Photon rate vs cos(Q)
TH2F * fMuonAnglesMi
Muon rate vs angle, middle momenta.
TH2F * fPhotonAnglesMi
Photon rate vs angle, middle momenta.
TH2F * fPhotonAnglesLo
Photon rate vs angle, low momenta.
TH2F * fPhotonAnglesHi
Photon rate vs angle, high momenta.
TH1F * fMuonsPerSample
number of muons in the sampled time window
TH1F * fMuonEnergy
Muon energy (GeV)
void evgen::CosmicsGen::beginRun ( art::Run run)
overrideprivatevirtual

Reimplemented from art::EDProducer.

Definition at line 142 of file CosmicsGen_module.cc.

143  {
145  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
146  }
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
LArSoft geometry interface.
Definition: ChannelGeo.h:16
void evgen::CosmicsGen::produce ( art::Event evt)
overrideprivatevirtual

Implements art::EDProducer.

Definition at line 149 of file CosmicsGen_module.cc.

150  {
151  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
152 
153  // fill some histograms about this event
155 
156  int nCrossCryostat = 0;
157 
158  simb::MCTruth truth;
159 
160  while(nCrossCryostat < 1){
161 
162  simb::MCTruth pretruth;
164  fCRYHelp.Sample(pretruth,
165  geom->SurfaceY(),
166  geom->DetLength(),
167  0);
168 
169  int numPhotons = 0;
170  int numElectrons = 0;
171  int numMuons = 0;
172  int allPhotons = 0;
173  int allElectrons = 0;
174  int allMuons = 0;
175 
176  // loop over particles in the truth object
177  for(int i = 0; i < pretruth.NParticles(); ++i){
178  simb::MCParticle particle = pretruth.GetParticle(i);
179  const TLorentzVector& v4 = particle.Position();
180  const TLorentzVector& p4 = particle.Momentum();
181  double x0[3] = {v4.X(), v4.Y(), v4.Z() };
182  double dx[3] = {p4.Px(), p4.Py(), p4.Pz()};
183 
184  if (std::abs(particle.PdgCode())==13) ++allMuons;
185  else if (std::abs(particle.PdgCode())==22) ++allPhotons;
186  else if (std::abs(particle.PdgCode())==11) ++allElectrons;
187 
188  TH1F* hCosQ = 0;
189  TH2F* hAngles = 0;
190  TH2F* hAnglesLo = 0;
191  TH2F* hAnglesMi = 0;
192  TH2F* hAnglesHi = 0;
193  TH1F* hEnergy = 0;
194  if (std::abs(particle.PdgCode())==13) {
195  hCosQ = fMuonCosQ;
196  hAngles = fMuonAngles;
197  hAnglesLo = fMuonAnglesLo;
198  hAnglesMi = fMuonAnglesMi;
199  hAnglesHi = fMuonAnglesHi;
200  hEnergy = fMuonEnergy;
201  }
202  else if (std::abs(particle.PdgCode())==22) {
203  hCosQ = fPhotonCosQ;
204  hAngles = fPhotonAngles;
205  hAnglesLo = fPhotonAnglesLo;
206  hAnglesMi = fPhotonAnglesMi;
207  hAnglesHi = fPhotonAnglesHi;
208  hEnergy = fPhotonEnergy;
209  }
210  else if (std::abs(particle.PdgCode())==11) {
211  hCosQ = fElectronCosQ;
212  hAngles = fElectronAngles;
213  hAnglesLo = fElectronAnglesLo;
214  hAnglesMi = fElectronAnglesMi;
215  hAnglesHi = fElectronAnglesHi;
216  hEnergy = fElectronEnergy;
217  }
218 
219  // now check if the particle goes through any cryostat in the detector
220  // if so, add it to the truth object.
221  for(unsigned int c = 0; c < geom->Ncryostats(); ++c){
222 
223  double bounds[6] = {0.};
224  geom->CryostatBoundaries(bounds, c);
225 
226  //add a buffer box around the cryostat bounds to increase the acceptance
227  //(geometrically) at the CRY level to make up for particles we will loose
228  //due to multiple scattering effects that pitch in during GEANT4 tracking
229  //By default, the buffer box has zero size
230  for (unsigned int cb=0; cb<6; cb++)
231  bounds[cb] = bounds[cb]+fbuffbox[cb];
232 
233  //calculate the intersection point with each cryostat surface
234  bool intersects_cryo = false;
235  for (int bnd=0; bnd!=6; ++bnd) {
236  if (bnd<2) {
237  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])};
238  if ( p2[1] >= bounds[2] && p2[1] <= bounds[3] &&
239  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
240  intersects_cryo = true;
241  break;
242  }
243  }
244  else if (bnd>=2 && bnd<4) {
245  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])};
246  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
247  p2[2] >= bounds[4] && p2[2] <= bounds[5] ) {
248  intersects_cryo = true;
249  break;
250  }
251  }
252  else if (bnd>=4) {
253  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]};
254  if ( p2[0] >= bounds[0] && p2[0] <= bounds[1] &&
255  p2[1] >= bounds[2] && p2[1] <= bounds[3] ) {
256  intersects_cryo = true;
257  break;
258  }
259  }
260  }
261 
262 
263  if (intersects_cryo) {
264  truth.Add(particle);
265 
266  if (std::abs(particle.PdgCode())==13) ++numMuons;
267  else if (std::abs(particle.PdgCode())==22) ++numPhotons;
268  else if (std::abs(particle.PdgCode())==11) ++numElectrons;
269 
270  //The following code no longer works now that we require intersection with the cryostat boundary
271  //For example, the particle could intersect this cryostat but miss its TPC, but intersect a TPC
272  //in another cryostat
273  /*try{
274  unsigned int tpc = 0;
275  unsigned int cstat = 0;
276  geom->PositionToTPC(x2, tpc, cstat);
277  if (std::abs(particle.PdgCode())==13) ++tpcMuons;
278  else if (std::abs(particle.PdgCode())==22) ++tpcPhotons;
279  else if (std::abs(particle.PdgCode())==11) ++tpcElectrons;
280  }
281  catch(cet::exception &e){
282  MF_LOG_DEBUG("CosmicsGen") << "current particle does not go through any tpc";
283  }*///
284 
285  if (hCosQ!=0) {
286  double cosq = -p4.Py()/p4.P();
287  double phi = std::atan2(p4.Pz(),p4.Px());
288  phi *= 180/M_PI;
289  hCosQ->Fill(cosq);
290  hAngles->Fill(phi,cosq);
291  if (p4.E()<1.0) hAnglesLo->Fill(phi,cosq);
292  else if (p4.E()<10.0) hAnglesMi->Fill(phi,cosq);
293  else hAnglesHi->Fill(phi,cosq);
294  hEnergy->Fill(p4.E());
295  }//end if there is a cos(theta) histogram
296  break; //leave loop over cryostats to avoid adding particle multiple times
297  }// end if particle goes into a cryostat
298  }// end loop over cryostats in the detector
299 
300  }// loop on particles
301 
302  nCrossCryostat = truth.NParticles();
303 
304  fPhotonsPerSample ->Fill(allPhotons);
305  fElectronsPerSample->Fill(allElectrons);
306  fMuonsPerSample ->Fill(allMuons);
307 
308  fPhotonsInCStat ->Fill(numPhotons);
309  fElectronsInCStat->Fill(numElectrons);
310  fMuonsInCStat ->Fill(numMuons);
311 
312  /*fPhotonsInTPC ->Fill(tpcPhotons);
313  fElectronsInTPC->Fill(tpcElectrons);
314  fMuonsInTPC ->Fill(tpcMuons);*/
315  }
316 
317  truthcol->push_back(truth);
318  evt.put(std::move(truthcol));
319  }// end produce
TH2F * fElectronAnglesHi
Electron rate vs angle, high momenta.
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:219
TH2F * fPhotonAngles
Photon rate vs angle.
int PdgCode() const
Definition: MCParticle.h:212
TH2F * fMuonAnglesHi
Muon rate vs angle, high momenta.
TH1F * fElectronEnergy
Electron energy (GeV)
TH2F * fElectronAnglesLo
Electron rate vs angle, low momenta.
TH2F * fElectronAngles
Electron rate vs angle.
TH1F * fElectronCosQ
Electron rate vs cos(Q)
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
TH2F * fMuonAngles
Muon rate vs angle.
TH2F * fMuonAnglesLo
Muon rate vs angle, low momenta.
TH1F * fPhotonEnergy
Photon energy (GeV)
int NParticles() const
Definition: MCTruth.h:75
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
TH2F * fElectronAnglesMi
Electron rate vs angle, middle momenta.
TH1F * fMuonCosQ
Muon rate vs cos(Q)
TH1F * fElectronsPerSample
number of electrons in the sampled time window
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
TH1F * fPhotonsPerSample
number of photons in the sampled time window
T abs(T value)
double Sample(simb::MCTruth &mctruth, double const &surfaceY, double const &detectorLength, double *w, double rantime=0)
Definition: CRYHelper.cxx:97
def move(depos, offset)
Definition: depos.py:107
TH1F * fPhotonCosQ
Photon rate vs cos(Q)
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
geo::Length_t SurfaceY() const
The position of the detector respect to earth surface.
evgb::CRYHelper fCRYHelp
CRY generator object.
BoundingBox bounds(int x, int y, int w, int h)
Definition: main.cpp:37
TH2F * fMuonAnglesMi
Muon rate vs angle, middle momenta.
std::vector< double > fbuffbox
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:76
TH2F * fPhotonAnglesMi
Photon rate vs angle, middle momenta.
TH2F * fPhotonAnglesLo
Photon rate vs angle, low momenta.
#define M_PI
Definition: includeROOT.h:54
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
TH2F * fPhotonAnglesHi
Photon rate vs angle, high momenta.
Event generator information.
Definition: MCTruth.h:32
TH1F * fMuonsPerSample
number of muons in the sampled time window
TH1F * fMuonEnergy
Muon energy (GeV)
Cosmic rays.
Definition: MCTruth.h:24

Member Data Documentation

std::vector<double> evgen::CosmicsGen::fbuffbox
private

Definition at line 46 of file CosmicsGen_module.cc.

evgb::CRYHelper evgen::CosmicsGen::fCRYHelp
private

CRY generator object.

Definition at line 84 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fElectronAngles
private

Electron rate vs angle.

Definition at line 60 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fElectronAnglesHi
private

Electron rate vs angle, high momenta.

Definition at line 63 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fElectronAnglesLo
private

Electron rate vs angle, low momenta.

Definition at line 61 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fElectronAnglesMi
private

Electron rate vs angle, middle momenta.

Definition at line 62 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fElectronCosQ
private

Electron rate vs cos(Q)

Definition at line 64 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fElectronEnergy
private

Electron energy (GeV)

Definition at line 65 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fElectronsInCStat
private

number of electrons in the cryostat during the sampled time window

Definition at line 67 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fElectronsInTPC
private

number of electrons in the tpc during the sampled time window

Definition at line 69 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fElectronsPerSample
private

number of electrons in the sampled time window

Definition at line 66 of file CosmicsGen_module.cc.

CLHEP::HepRandomEngine& evgen::CosmicsGen::fEngine
private

art-managed random-number engine

Definition at line 83 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fMuonAngles
private

Muon rate vs angle.

Definition at line 72 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fMuonAnglesHi
private

Muon rate vs angle, high momenta.

Definition at line 75 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fMuonAnglesLo
private

Muon rate vs angle, low momenta.

Definition at line 73 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fMuonAnglesMi
private

Muon rate vs angle, middle momenta.

Definition at line 74 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fMuonCosQ
private

Muon rate vs cos(Q)

Definition at line 76 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fMuonEnergy
private

Muon energy (GeV)

Definition at line 77 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fMuonsInCStat
private

number of muons in the cryostat during the sampled time window

Definition at line 79 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fMuonsInTPC
private

number of muons in the tpc during the sampled time window

Definition at line 81 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fMuonsPerSample
private

number of muons in the sampled time window

Definition at line 78 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fPhotonAngles
private

Photon rate vs angle.

Definition at line 48 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fPhotonAnglesHi
private

Photon rate vs angle, high momenta.

Definition at line 51 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fPhotonAnglesLo
private

Photon rate vs angle, low momenta.

Definition at line 49 of file CosmicsGen_module.cc.

TH2F* evgen::CosmicsGen::fPhotonAnglesMi
private

Photon rate vs angle, middle momenta.

Definition at line 50 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fPhotonCosQ
private

Photon rate vs cos(Q)

Definition at line 52 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fPhotonEnergy
private

Photon energy (GeV)

Definition at line 53 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fPhotonsInCStat
private

number of photons in the cryostat during the sampled time window

Definition at line 55 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fPhotonsInTPC
private

number of photons in the tpc during the sampled time window

Definition at line 57 of file CosmicsGen_module.cc.

TH1F* evgen::CosmicsGen::fPhotonsPerSample
private

number of photons in the sampled time window

Definition at line 54 of file CosmicsGen_module.cc.


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