All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SignalShapingServiceSBND.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////
2 ///
3 /// \file SignalShapingServiceSBND.h
4 ///
5 /// \brief Service to provide microboone-specific signal shaping for
6 /// simulation (convolution) and reconstruction (deconvolution).
7 ///
8 /// \author H. Greenlee
9 ///
10 /// This service inherits from SignalShaping and supplies
11 /// microboone-specific configuration. It is intended that SimWire and
12 /// CalWire modules will access this service.
13 ///
14 /// FCL parameters:
15 ///
16 /// FieldBins - Number of bins of field response.
17 /// Col3DCorrection - 3D path length correction for collection plane.
18 /// Ind3DCorrection - 3D path length correction for induction plane.
19 /// ColFieldRespAmp - Collection field response amplitude.
20 /// IndFieldRespAmp - Induction field response amplitude.
21 /// ShapeTimeConst - Time constants for exponential shaping.
22 /// ColFilter - Root parameterized collection plane filter function.
23 /// ColFilterParams - Collection filter function parameters.
24 /// IndFilter - Root parameterized induction plane filter function.
25 /// IndFilterParams - Induction filter function parameters.
26 ///
27 ////////////////////////////////////////////////////////////////////////
28 
29 #ifndef SIGNALSHAPINGSERVICELARIAT_H
30 #define SIGNALSHAPINGSERVICELARIAT_H
31 
32 #include <vector>
33 
34 #include "fhiclcpp/ParameterSet.h"
35 #include "art/Framework/Services/Registry/ActivityRegistry.h"
36 #include "art/Framework/Services/Registry/ServiceDeclarationMacros.h"
41 namespace detinfo { class DetectorClocksData; }
42 
43 #include "TF1.h"
44 #include "TH1D.h"
45 
46 
47 using DoubleVec = std::vector<double>;
48 
49 namespace util {
51  public:
52 
53  // Constructor, destructor.
54 
55  SignalShapingServiceSBND(const fhicl::ParameterSet& pset,
56  art::ActivityRegistry& reg);
58 
59  // Update configuration parameters.
60 
61  void reconfigure(const fhicl::ParameterSet& pset);
62 
63  std::vector<DoubleVec> GetNoiseFactVec() {return fNoiseFactVec;};
64  double GetASICGain(unsigned int const channel) const;
65  //double GetShapingTime(unsigned int const channel) const;
66 
67  double GetRawNoise(unsigned int const channel) const;
68  double GetDeconNoise(unsigned int const channel) const;
69 
70  // Accessors.
71 
72  const util::SignalShaping& SignalShaping(unsigned int channel) const;
73 
75  unsigned int const channel) const;
76 
77  // Do convolution calcution (for simulation).
78 
79  template <class T> void Convolute(detinfo::DetectorClocksData const& clockData,
80  unsigned int channel, std::vector<T>& func) const;
81 
82  // Do deconvolution calcution (for reconstruction).
83 
84  template <class T> void Deconvolute(detinfo::DetectorClocksData const& clockData,
85  unsigned int channel, std::vector<T>& func) const;
86 
87  double GetDeconNorm(){return fDeconNorm;};
88 
89  private:
90 
91  // Private configuration methods.
92 
93  // Post-constructor initialization.
94 
95  void init() const{const_cast<SignalShapingServiceSBND*>(this)->init();}
96  void init();
97 
98 
99  // Calculate response functions.
100  // Copied from SimWireSBND.
101 
102  void SetFieldResponse();
103  void SetElectResponse(double shapingtime, double gain);
104 
105  // Calculate filter functions.
106 
107 
108  void SetFilters();
109 
110  // Calculate view corresponding to channel
111  geo::View_t GetView(unsigned int chan) const;
112 
113  // Attributes.
114 
115  bool fInit; ///< Initialization flag.
116 
117  void SetResponseSampling();
118 
119  // Fcl parameters.
120  double fDeconNorm;
121  double fADCPerPCAtLowestASICGain; ///Pulse amplitude gain for a 1 pc charge impulse after convoluting it with field and electronics response with the lowest ASIC gain setting of 4.7 mV/fC
122 
123  std::vector<DoubleVec> fNoiseFactVec;
124  std::vector<double> fASICGainInMVPerFC;
125  std::vector<double> fShapeTimeConst; ///< time constants for exponential shaping
126 
127  std::vector<double> fFieldResponseTOffset; ///<Time offset for field response in ns
128  double fInputFieldRespSamplingPeriod; ///< Sampling Period in the input field response
129 
130  int fNFieldBins; ///< number of bins for field response
131  double fCol3DCorrection; ///< correction factor to account for 3D path of
132  ///< electrons thru wires
133  double fInd3DCorrection; ///< correction factor to account for 3D path of
134  ///< electrons thru wires
135  double fColFieldRespAmp; ///< amplitude of response to field
136  double fIndUFieldRespAmp; ///< amplitude of response to field in U plane
137  double fIndVFieldRespAmp; ///< amplitude of response to field in V plane
138  TF1* fColFilterFunc; ///< Parameterized collection filter function.
139  TF1* fIndUFilterFunc; ///< Parameterized induction filter function for U plane.
140  TF1* fIndVFilterFunc; ///< Parameterized induction filter function for V plane
141 
142 
143  bool fUseFunctionFieldShape; ///< Flag that allows to use a parameterized field response instead of the hardcoded version
144  bool fUseSimpleFieldShape; ///< Flag that turns on new field response shapes
145  bool fUseHistogramFieldShape; ///< Flag that truns on field response shaped from histograms
146  bool fGetFilterFromHisto; ///< Flag that allows to use a filter function from a histogram instead of the functional dependency
147  TF1* fColFieldFunc; ///< Parameterized collection field shape function.
148  TF1* fIndUFieldFunc; ///< Parameterized induction field shape function for U plane.
149  TF1* fIndVFieldFunc; ///< Parameterized induction field shape function for V plane.
150 
151  TH1 const* fFieldResponseHist[3]; ///< Histogram used to hold the field response, hardcoded for the time being
152  TH1 const* fFilterHist[3]; ///< Histogram used to hold the collection filter, hardcoded for the time being
153 
154  // Following attributes hold the convolution and deconvolution kernels
155 
159 
160  // Field response.
161 
162  std::vector<double> fIndUFieldResponse;
163  std::vector<double> fIndVFieldResponse;
164  std::vector<double> fColFieldResponse;
165 
166  // Electronics response.
167 
168  std::vector<double> fElectResponse;
169 
170  // Filters.
171 
172  std::vector<TComplex> fIndUFilter;
173  std::vector<TComplex> fIndVFilter;
174  std::vector<TComplex> fColFilter;
175  };
176 }
177 //----------------------------------------------------------------------
178 // Do convolution.
179 template <class T> inline void util::SignalShapingServiceSBND::Convolute(detinfo::DetectorClocksData const& clockData,
180  unsigned int channel, std::vector<T>& func) const
181 {
182  SignalShaping(channel).Convolute(func);
183 
184  //negative number;
185  int time_offset = FieldResponseTOffset(clockData, channel);
186 
187  std::vector<T> temp;
188  if (time_offset <= 0) {
189  temp.assign(func.begin(),func.begin()-time_offset);
190  func.erase(func.begin(),func.begin()-time_offset);
191  func.insert(func.end(),temp.begin(),temp.end());
192  }else{
193  temp.assign(func.end()-time_offset,func.end());
194  func.erase(func.end()-time_offset,func.end());
195  func.insert(func.begin(),temp.begin(),temp.end());
196  }
197 }
198 
199 
200 //----------------------------------------------------------------------
201 // Do deconvolution.
202 template <class T> inline void util::SignalShapingServiceSBND::Deconvolute(detinfo::DetectorClocksData const& clockData,
203  unsigned int channel, std::vector<T>& func) const
204 {
205  SignalShaping(channel).Deconvolute(func);
206 
207  //negative number;
208  int time_offset = FieldResponseTOffset(clockData, channel);
209 
210  std::vector<T> temp;
211  if (time_offset <= 0) {
212  temp.assign(func.end()+time_offset,func.end());
213  func.erase(func.end()+time_offset,func.end());
214  func.insert(func.begin(),temp.begin(),temp.end());
215  }else{
216  temp.assign(func.begin(),func.begin()+time_offset);
217  func.erase(func.begin(),func.begin()+time_offset);
218  func.insert(func.end(),temp.begin(),temp.end());
219  }
220 
221 }
222 
223 DECLARE_ART_SERVICE(util::SignalShapingServiceSBND, LEGACY)
224 #endif
double GetDeconNoise(unsigned int const channel) const
void Convolute(detinfo::DetectorClocksData const &clockData, unsigned int channel, std::vector< T > &func) const
std::vector< double > DoubleVec
void Deconvolute(std::vector< T > &func) const
TH1 const * fFilterHist[3]
Histogram used to hold the collection filter, hardcoded for the time being.
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
std::vector< double > fShapeTimeConst
time constants for exponential shaping
TF1 * fColFilterFunc
Parameterized collection filter function.
double fIndUFieldRespAmp
amplitude of response to field in U plane
double GetRawNoise(unsigned int const channel) const
void SetElectResponse(double shapingtime, double gain)
double fInputFieldRespSamplingPeriod
Sampling Period in the input field response.
TF1 * fIndVFieldFunc
Parameterized induction field shape function for V plane.
int fNFieldBins
number of bins for field response
TF1 * fColFieldFunc
Parameterized collection field shape function.
void Convolute(std::vector< T > &func) const
bool fUseSimpleFieldShape
Flag that turns on new field response shapes.
Generic class for shaping signals on wires.
std::vector< double > fFieldResponseTOffset
Time offset for field response in ns.
bool fUseFunctionFieldShape
Flag that allows to use a parameterized field response instead of the hardcoded version.
const util::SignalShaping & SignalShaping(unsigned int channel) const
TF1 * fIndUFieldFunc
Parameterized induction field shape function for U plane.
bool fUseHistogramFieldShape
Flag that truns on field response shaped from histograms.
double GetASICGain(unsigned int const channel) const
TF1 * fIndUFilterFunc
Parameterized induction filter function for U plane.
std::vector< DoubleVec > fNoiseFactVec
Pulse amplitude gain for a 1 pc charge impulse after convoluting it with field and electronics respon...
TH1 const * fFieldResponseHist[3]
Histogram used to hold the field response, hardcoded for the time being.
double fColFieldRespAmp
amplitude of response to field
SignalShapingServiceSBND(const fhicl::ParameterSet &pset, art::ActivityRegistry &reg)
double fIndVFieldRespAmp
amplitude of response to field in V plane
Encapsulate the construction of a single detector plane.
Contains all timing reference information for the detector.
bool fGetFilterFromHisto
Flag that allows to use a filter function from a histogram instead of the functional dependency...
int FieldResponseTOffset(detinfo::DetectorClocksData const &clockData, unsigned int const channel) const
void Deconvolute(detinfo::DetectorClocksData const &clockData, unsigned int channel, std::vector< T > &func) const
std::vector< DoubleVec > GetNoiseFactVec()
geo::View_t GetView(unsigned int chan) const
TF1 * fIndVFilterFunc
Parameterized induction filter function for V plane.
art framework interface to geometry description
Encapsulate the construction of a single detector plane.
void reconfigure(const fhicl::ParameterSet &pset)