All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Attributes | List of all members
SBNDNoiseServiceFromHist Class Reference

#include <SBNDNoiseServiceFromHist.h>

Inheritance diagram for SBNDNoiseServiceFromHist:
ChannelNoiseService

Public Member Functions

 SBNDNoiseServiceFromHist (fhicl::ParameterSet const &pset)
 
 SBNDNoiseServiceFromHist (fhicl::ParameterSet const &pset, art::ActivityRegistry &)
 
 ~SBNDNoiseServiceFromHist ()
 
int addNoise (Channel chan, AdcSignalVector &sigs) const
 
std::ostream & print (std::ostream &out=std::cout, std::string prefix="") const
 
- Public Member Functions inherited from ChannelNoiseService
virtual ~ChannelNoiseService ()=default
 
virtual int addNoise (detinfo::DetectorClocksData const &, Channel chan, AdcSignalVector &sigs) const =0
 
virtual void generateNoise (detinfo::DetectorClocksData const &)
 
virtual void InitialiseProducerDeps (art::EDProducer *EDProdPointer, fhicl::ParameterSet const &pset)
 

Private Attributes

unsigned int fNoiseArrayPoints
 number of points in randomly generated noise array More...
 
int fRandomSeed
 Seed for random number service. If absent or zero, use SeedSvc. More...
 
int fLogLevel
 Log message level: 0=quiet, 1=init only, 2+=every event. More...
 
CLHEP::HepRandomEngine * m_pran
 

Additional Inherited Members

- Public Types inherited from ChannelNoiseService
typedef unsigned int Channel
 

Detailed Description

Definition at line 51 of file SBNDNoiseServiceFromHist.h.

Constructor & Destructor Documentation

SBNDNoiseServiceFromHist::SBNDNoiseServiceFromHist ( fhicl::ParameterSet const &  pset)

Definition at line 20 of file SBNDNoiseServiceFromHist_service.cc.

21  : fRandomSeed(0), fLogLevel(1), m_pran(nullptr), fNoiseEngine(nullptr)
22 {
23  const string myname = "SBNDNoiseServiceFromHist::ctor: ";
24  fNoiseArrayPoints = pset.get<unsigned int>("NoiseArrayPoints");
25  bool haveSeed = pset.get_if_present<int>("RandomSeed", fRandomSeed);
26 
27 
28  fShapingTimeOrder = { {0.5, 0 }, {1.0, 1}, {2.0, 2}, {3.0, 3} };
29  fNoiseWidth = pset.get< double >("NoiseWidth");
30  fNoiseRand = pset.get< double >("NoiseRand");
31  fLowCutoff = pset.get< double >("LowCutoff");
32 
33  //Getting noise histo
34  fNoiseHistoName = p.get< std::string >("NoiseHistoName");
35 
36  cet::search_path sp("FW_SEARCH_PATH");
37  sp.find_file(p.get<std::string>("NoiseFileFname"), fNoiseFileFname);
38 
39  TFile in(fNoiseFileFname.c_str(), "READ");
40  if (!in.IsOpen()) {
41  throw art::Exception(art::errors::FileOpenError)
42  << "Could not find Root file '" << fNoiseHistoName
43  << "' for noise histogram\n";
44  }
45  fNoiseHist = (TH1D*) in.Get(fNoiseHistoName.c_str());
46 
47  if (!fNoiseHist) {
48  throw art::Exception(art::errors::NotFound)
49  << "Could not find noise histogram '" << fNoiseHistoName
50  << "' in Root file '" << in.GetPath() << "'\n";
51  }
52  // release the histogram from its file; now we own it
53  fNoiseHist->SetDirectory(nullptr);
54  in.Close();
55 
56 
57  if ( fRandomSeed == 0 ) haveSeed = false;
58  pset.get_if_present<int>("LogLevel", fLogLevel);
59  int seed = fRandomSeed;
60 
61  string rname = "SBNDNoiseServiceFromHist";
62  if ( haveSeed ) {
63  if ( fLogLevel > 0 ) cout << myname << "WARNING: Using hardwired seed." << endl;
64  m_pran = new HepJamesRandom(seed);
65  } else {
66  if ( fLogLevel > 0 ) cout << myname << "Using NuRandomService." << endl;
67  art::ServiceHandle<NuRandomService> seedSvc;
68  m_pran = new HepJamesRandom;
69  if ( fLogLevel > 0 ) cout << myname << " Initial seed: " << m_pran->getSeed() << endl;
70  seedSvc->registerEngine(NuRandomService::CLHEPengineSeeder(m_pran), rname);
71  }
72  if ( fLogLevel > 0 ) cout << myname << " Registered seed: " << m_pran->getSeed() << endl;
73  generateNoise();
74  if ( fLogLevel > 1 ) print() << endl;
75 }
CLHEP::HepRandomEngine * m_pran
pdgs p
Definition: selectors.fcl:22
int fRandomSeed
Seed for random number service. If absent or zero, use SeedSvc.
std::ostream & print(std::ostream &out=std::cout, std::string prefix="") const
int fLogLevel
Log message level: 0=quiet, 1=init only, 2+=every event.
unsigned int seed
if &&[-z"$BASH_VERSION"] then echo Attempting to switch to bash bash shellSwitch exit fi &&["$1"= 'shellSwitch'] shift declare a IncludeDirectives for Dir in
virtual void generateNoise(detinfo::DetectorClocksData const &)
unsigned int fNoiseArrayPoints
number of points in randomly generated noise array
BEGIN_PROLOG could also be cout
SBNDNoiseServiceFromHist::SBNDNoiseServiceFromHist ( fhicl::ParameterSet const &  pset,
art::ActivityRegistry &   
)

Definition at line 80 of file SBNDNoiseServiceFromHist_service.cc.

81 : SBNDNoiseServiceFromHist(pset) { }
SBNDNoiseServiceFromHist(fhicl::ParameterSet const &pset)
SBNDNoiseServiceFromHist::~SBNDNoiseServiceFromHist ( )

Definition at line 85 of file SBNDNoiseServiceFromHist_service.cc.

85  {
86  const string myname = "SBNDNoiseServiceFromHist::dtor: ";
87  if ( fLogLevel > 0 ) {
88  cout << myname << "Deleting random engine with seed " << m_pran->getSeed() << endl;
89  }
90  delete m_pran;
91 }
CLHEP::HepRandomEngine * m_pran
int fLogLevel
Log message level: 0=quiet, 1=init only, 2+=every event.
BEGIN_PROLOG could also be cout

Member Function Documentation

int SBNDNoiseServiceFromHist::addNoise ( Channel  chan,
AdcSignalVector sigs 
) const

Definition at line 95 of file SBNDNoiseServiceFromHist_service.cc.

95  {
96 
97  //Get services.
98  art::ServiceHandle<geo::Geometry> geo;
99  art::ServiceHandle<util::SignalShapingServiceSBND> sss;
100  art::ServiceHandle<util::LArFFT> fFFT;
101 
102  //Generate Noise:
103  size_t view = (size_t)geo->View(chan);
104 
105  double noise_factor;
106  auto tempNoiseVec = sss->GetNoiseFactVec();
107  double shapingTime = 2.0; //sss->GetShapingTime(chan);
108  double asicGain = sss->GetASICGain(chan);
109 
110  size_t fNTicks = fFFT->FFTSize();
111 
112 
113  if (fShapingTimeOrder.find( shapingTime ) != fShapingTimeOrder.end() ) {
114  noise_factor = tempNoiseVec[view].at( fShapingTimeOrder.find( shapingTime )->second );
115  noise_factor *= asicGain/4.7;
116  }
117  else {
118  throw cet::exception("SBNDNoiseServiceFromHist_service.cc")
119  << "\033[93m"
120  << "Shaping Time recieved from signalshapingservices_sbnd.fcl is not one of the allowed values"
121  << std::endl
122  << "Allowed values: 0.5, 1.0, 2.0, 3.0 us"
123  << "\033[00m"
124  << std::endl;
125  }
126 
127 
128  CLHEP::RandFlat flat(*fNoiseEngine, -1, 1);
129 
130 
131  if (sigs.size() != fNTicks)
132  throw cet::exception("SBNDNoiseServiceFromHist_service.cc")
133  << "\033[93m"
134  << "Frequency noise vector length must match fNTicks (FFT size)"
135  << " ... " << sigs.size() << " != " << fNTicks
136  << "\033[00m"
137  << std::endl;
138 
139  // noise in frequency space
140  std::vector<TComplex> noiseFrequency(fNTicks / 2 + 1, 0.);
141 
142  double pval = 0.;
143  double lofilter = 0.;
144  double phase = 0.;
145  double rnd[2] = {0.};
146 
147  // width of frequencyBin in kHz
148  double binWidth = 1.0 / (fNTicks * fSampleRate * 1.0e-6);
149 
150  for (size_t i = 0; i < fNTicks / 2 + 1; ++i) {
151  // exponential noise spectrum
152  flat.fireArray(2, rnd, 0, 1);
153  pval = fNoiseHist->GetBinContent(i) * ((1 - fNoiseRand) + 2 * fNoiseRand * rnd[0]) * noise_factor;
154  phase = rnd[1] * 2.*TMath::Pi();
155  TComplex tc(pval * cos(phase), pval * sin(phase));
156  noiseFrequency.at(i) += tc;
157  }
158 
159 
160  // inverse FFT MCSignal
161  fFFT->DoInvFFT(noiseFrequency, sigs);
162 
163  noiseFrequency.clear();
164 
165  // multiply each noise value by fNTicks as the InvFFT
166  // divides each bin by fNTicks assuming that a forward FFT
167  // has already been done.
168  for (unsigned int i = 0; i < sigs.size(); ++i) {
169  sigs.at(i) *= 1.*fNTicks;
170  }
171 
172  return 0;
173 }
ostream & SBNDNoiseServiceFromHist::print ( std::ostream &  out = std::cout,
std::string  prefix = "" 
) const
virtual

Implements ChannelNoiseService.

Definition at line 178 of file SBNDNoiseServiceFromHist_service.cc.

178  {
179  out << prefix << "SBNDNoiseServiceFromHist: " << endl;
180 
181  out << prefix << " LogLevel: " << fLogLevel << endl;
182  out << prefix << " RandomSeed: " << fRandomSeed << endl;
183  out << prefix << " NoiseArrayPoints: " << fNoiseArrayPoints << endl;
184 
185  return out;
186 }
int fRandomSeed
Seed for random number service. If absent or zero, use SeedSvc.
int fLogLevel
Log message level: 0=quiet, 1=init only, 2+=every event.
unsigned int fNoiseArrayPoints
number of points in randomly generated noise array

Member Data Documentation

int SBNDNoiseServiceFromHist::fLogLevel
private

Log message level: 0=quiet, 1=init only, 2+=every event.

Definition at line 75 of file SBNDNoiseServiceFromHist.h.

unsigned int SBNDNoiseServiceFromHist::fNoiseArrayPoints
private

number of points in randomly generated noise array

Definition at line 73 of file SBNDNoiseServiceFromHist.h.

int SBNDNoiseServiceFromHist::fRandomSeed
private

Seed for random number service. If absent or zero, use SeedSvc.

Definition at line 74 of file SBNDNoiseServiceFromHist.h.

CLHEP::HepRandomEngine* SBNDNoiseServiceFromHist::m_pran
private

Definition at line 77 of file SBNDNoiseServiceFromHist.h.


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