All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SBNDThermalNoiseServiceInTime_service.cc
Go to the documentation of this file.
1 // SBNDThermalNoiseServiceInTime.cxx
2 // Andrew Scarff
3 // July 2019
4 // Based upon SPhaseChannelNoiseService.cxx developed by Jingbo Wang for ProtoDUNE.
5 
7 #include "art/Framework/Services/Registry/ServiceDefinitionMacros.h"
8 
9 using std::cout;
10 using std::ostream;
11 using std::endl;
12 using std::string;
13 using std::ostringstream;
14 using rndm::NuRandomService;
15 using CLHEP::HepJamesRandom;
16 using CLHEP::HepRandomEngine;
17 
18 //**********************************************************************
19 
21 SBNDThermalNoiseServiceInTime(fhicl::ParameterSet const& pset): fRandomSeed(0), fLogLevel(1),
22  m_pran(nullptr), fNoiseEngine(nullptr)
23  {
24 
25  const string myname = "SBNDThermalNoiseServiceInTime::ctor: ";
26  fNoiseArrayPoints = pset.get<unsigned int>("NoiseArrayPoints",1);
27  bool haveSeed = pset.get_if_present<int>("RandomSeed", fRandomSeed);
28 
29  fShapingTimeOrder = { {0.5, 0 }, {1.0, 1}, {2.0, 2}, {3.0, 3} };
30 
31  if ( fRandomSeed == 0 ) haveSeed = false;
32  pset.get_if_present<int>("LogLevel", fLogLevel);
33  int seed = fRandomSeed;
34 
35  string rname = "SBNDThermalNoiseServiceInTime";
36  if ( haveSeed ) {
37  if ( fLogLevel > 0 ) cout << myname << "WARNING: Using hardwired seed." << endl;
38  m_pran = new HepJamesRandom(seed);
39  } else {
40  if ( fLogLevel > 0 ) cout << myname << "Using NuRandomService." << endl;
41  art::ServiceHandle<NuRandomService> seedSvc;
42  m_pran = new HepJamesRandom;
43  if ( fLogLevel > 0 ) cout << myname << " Initial seed: " << m_pran->getSeed() << endl;
44  seedSvc->registerEngine(NuRandomService::CLHEPengineSeeder(m_pran), rname);
45  }
46  if ( fLogLevel > 0 ) cout << myname << " Registered seed: " << m_pran->getSeed() << endl;
47  if ( fLogLevel > 1 ) print() << endl;
48 }
49 
50 //**********************************************************************
51 
53 SBNDThermalNoiseServiceInTime(fhicl::ParameterSet const& pset, art::ActivityRegistry&)
55 
56 //**********************************************************************
57 
59  const string myname = "SBNDThermalNoiseServiceInTime::dtor: ";
60  if ( fLogLevel > 0 ) {
61  cout << myname << "Deleting random engine with seed " << m_pran->getSeed() << endl;
62  }
63  delete m_pran;
64 }
65 
66 //**********************************************************************
67 
69  Channel chan, AdcSignalVector& sigs) const {
70 
71  //Get services.
72  art::ServiceHandle<geo::Geometry> geo;
73  art::ServiceHandle<util::SignalShapingServiceSBND> sss;
74 
75  //Generate Noise:
76  size_t view = (size_t)geo->View(chan);
77 
78  double noise_factor;
79  auto tempNoiseVec = sss->GetNoiseFactVec();
80  double shapingTime = 2.0; //sss->GetShapingTime(chan);
81  double asicGain = sss->GetASICGain(chan);
82 
83  if (fShapingTimeOrder.find( shapingTime ) != fShapingTimeOrder.end() ) {
84  noise_factor = tempNoiseVec[view].at( fShapingTimeOrder.find( shapingTime )->second );
85  noise_factor *= asicGain/4.7;
86  }
87  else {
88  throw cet::exception("SBNDThermalNoiseServiceInTime_service")
89  << "\033[93m"
90  << "Shaping Time recieved from signalshapingservices_sbnd.fcl is not one of the allowed values"
91  << std::endl
92  << "Allowed values: 0.5, 1.0, 2.0, 3.0 us"
93  << "\033[00m"
94  << std::endl;
95  }
96 
97  CLHEP::RandGaussQ rGauss(*fNoiseEngine, 0.0, noise_factor);
98 
99 
100  //In this case fNoiseFact is a value in ADC counts
101  //It is going to be the Noise RMS
102  //loop over all bins in "noise" vector
103  //and insert random noise value
104 
105  for (unsigned int i = 0; i < sigs.size(); i++){
106  sigs.at(i) = rGauss.fire();
107  }
108 
109  return 0;
110 }
111 
112 
113 //**********************************************************************
114 
115 ostream& SBNDThermalNoiseServiceInTime::print(ostream& out, string prefix) const {
116  out << prefix << "SBNDThermalNoiseServiceInTime: " << endl;
117 
118  out << prefix << " LogLevel: " << fLogLevel << endl;
119  out << prefix << " RandomSeed: " << fRandomSeed << endl;
120  out << prefix << " NoiseArrayPoints: " << fNoiseArrayPoints << endl;
121 
122  return out;
123 }
124 
125 
126 
127 //**********************************************************************
128 
129 DEFINE_ART_SERVICE_INTERFACE_IMPL(SBNDThermalNoiseServiceInTime, ChannelNoiseService)
130 
131 //**********************************************************************
int fLogLevel
Log message level: 0=quiet, 1=init only, 2+=every event.
std::ostream & print(std::ostream &out=std::cout, std::string prefix="") const override
SBNDThermalNoiseServiceInTime(fhicl::ParameterSet const &pset)
unsigned int fNoiseArrayPoints
number of points in randomly generated noise array
int addNoise(detinfo::DetectorClocksData const &clockData, Channel chan, AdcSignalVector &sigs) const override
unsigned int seed
std::vector< AdcSignal > AdcSignalVector
Contains all timing reference information for the detector.
int fRandomSeed
Seed for random number service. If absent or zero, use SeedSvc.
BEGIN_PROLOG could also be cout