Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CaloCalibration.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file CaloCalibration.cc
1 #include "CaloCalibration.h"
2 
3 #include "PROTOTYPE4_FEM.h"
4 #include "RawTower_Prototype4.h"
5 
6 #include <calobase/RawTower.h> // for RawTower
7 #include <calobase/RawTowerContainer.h>
8 #include <calobase/RawTowerDefs.h> // for keytype
9 
10 #include <phparameter/PHParameters.h> // for PHParameters
11 
13 #include <fun4all/SubsysReco.h> // for SubsysReco
14 
15 #include <phool/PHCompositeNode.h>
16 #include <phool/PHIODataNode.h> // for PHIODataNode
17 #include <phool/PHNode.h> // for PHNode
18 #include <phool/PHNodeIterator.h> // for PHNodeIterator
19 #include <phool/PHObject.h> // for PHObject
20 #include <phool/getClass.h>
21 
22 #include <boost/format.hpp>
23 
24 #include <cassert>
25 #include <cmath> // for NAN, isnan
26 #include <cstdlib> // for exit
27 #include <iostream>
28 #include <limits> // for numeric_limits
29 #include <map> // for map, _Rb_tree_iter...
30 #include <stdexcept> // for runtime_error
31 #include <string>
32 #include <utility> // for pair
33 #include <vector> // for vector
34 
35 using namespace std;
36 
37 //____________________________________
39  : SubsysReco(string("CaloCalibration_") + name)
40  , _calib_towers(nullptr)
41  , _raw_towers(nullptr)
42  , detector(name)
43  , _calib_tower_node_prefix("CALIB")
44  , _raw_tower_node_prefix("RAW")
45  , _calib_params(name)
46  , _fit_type(kPowerLawDoubleExpWithGlobalFitConstraint)
47 {
49 }
50 
51 //_____________________________________
53 {
54  CreateNodeTree(topNode);
55 
56  if (Verbosity())
57  {
58  std::cout << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
59  << " - print calibration parameters: " << endl;
61  }
62 
64 }
65 
66 //____________________________________
68 {
69  if (Verbosity())
70  {
71  std::cout << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
72  << "Process event entered" << std::endl;
73  }
74 
75  map<int, double> parameters_constraints;
77  _raw_towers->size() > 1)
78  {
79  if (Verbosity())
80  {
81  std::cout
82  << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
83  << "Extract global fit parameter for constraining individual fits"
84  << std::endl;
85  }
86 
87  // signal template
88 
89  vector<double> vec_signal_samples(PROTOTYPE4_FEM::NSAMPLES, 0);
90 
91  int count = 0;
92 
95  for (rtiter = begin_end.first; rtiter != begin_end.second; ++rtiter)
96  {
97  RawTower_Prototype4 *raw_tower =
98  dynamic_cast<RawTower_Prototype4 *>(rtiter->second);
99  assert(raw_tower);
100 
101  // bool signal_check_pass = true;
102  // for (int i = 0; i < RawTower_Prototype4::NSAMPLES; i++)
103  // {
104  // if (raw_tower->get_signal_samples(i) <= 10 or
105  // raw_tower->get_signal_samples(i) >= ((1 << 14) - 10))
106  // {
107  // signal_check_pass = false;
108  // break;
109  // }
110  // }
111 
112  // if (signal_check_pass)
113  // {
114  ++count;
115 
116  for (int i = 0; i < RawTower_Prototype4::NSAMPLES; i++)
117  {
118  if (raw_tower->get_signal_samples(i) <= 10 or
119  raw_tower->get_signal_samples(i) >= ((1 << 14) - 10))
120  vec_signal_samples[i] =
121  numeric_limits<double>::quiet_NaN(); // invalidate this sample
122  else
123  vec_signal_samples[i] += raw_tower->get_signal_samples(i);
124  }
125  // }
126  }
127 
128  if (count > 0)
129  {
130  for (int i = 0; i < RawTower_Prototype4::NSAMPLES; i++)
131  {
132  vec_signal_samples[i] /= count;
133  }
134 
135  double peak = NAN;
136  double peak_sample = NAN;
137  double pedstal = NAN;
138  map<int, double> parameters_io;
139 
140  PROTOTYPE4_FEM::SampleFit_PowerLawDoubleExp(vec_signal_samples, peak,
141  peak_sample, pedstal,
142  parameters_io, Verbosity());
143  // std::map<int, double> &parameters_io, //! IO for fullset of
144  // parameters. If a parameter exist and not an NAN, the fit parameter
145  // will be fixed to that value. The order of the parameters are
146  // ("Amplitude", "Sample Start", "Power", "Peak Time 1", "Pedestal",
147  // "Amplitude ratio", "Peak Time 2")
148 
149  parameters_constraints[1] = parameters_io[1];
150  parameters_constraints[2] = parameters_io[2];
151  parameters_constraints[3] = parameters_io[3];
152  parameters_constraints[5] = parameters_io[5];
153  parameters_constraints[6] = parameters_io[6];
154 
155  // //special constraint if Peak Time 1 == Peak Time 2
156  // if (abs(parameters_constraints[6] - parameters_constraints[3]) <
157  // 0.1)
158  // {
159  // const double average_peak_time = (parameters_constraints[6] +
160  // parameters_constraints[3]) / 2.;
161  //
162  // std::cout << Name() << "::" << detector << "::" <<
163  // __PRETTY_FUNCTION__
164  // << ": two shaping time are too close "
165  // << parameters_constraints[3] << " VS " <<
166  // parameters_constraints[6]
167  // << ". Use average peak time instead: " <<
168  // average_peak_time
169  // << std::endl;
170  //
171  // parameters_constraints[6] = average_peak_time;
172  // parameters_constraints[3] = average_peak_time;
173  // parameters_constraints[5] = 0;
174  // }
175  }
176  else
177  {
178  std::cout << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
179  << ": Failed to build signal template! Fit each channel "
180  "individually instead"
181  << std::endl;
182  }
183  }
184 
185  const double calib_const_scale =
186  _calib_params.get_double_param("calib_const_scale");
187  const bool use_chan_calibration =
188  _calib_params.get_int_param("use_chan_calibration") > 0;
189 
192  for (rtiter = begin_end.first; rtiter != begin_end.second; ++rtiter)
193  {
194  RawTowerDefs::keytype key = rtiter->first;
195  RawTower_Prototype4 *raw_tower =
196  dynamic_cast<RawTower_Prototype4 *>(rtiter->second);
197  assert(raw_tower);
198 
199  double calibration_const = calib_const_scale;
200 
201  if (use_chan_calibration)
202  {
203  // channel to channel calibration.
204  const int column = raw_tower->get_column();
205  const int row = raw_tower->get_row();
206 
207  assert(column >= 0);
208  assert(row >= 0);
209 
210  string calib_const_name(
211  (boost::format("calib_const_column%d_row%d") % column % row).str());
212 
213  calibration_const *= _calib_params.get_double_param(calib_const_name);
214  }
215 
216  vector<double> vec_signal_samples;
217  for (int i = 0; i < RawTower_Prototype4::NSAMPLES; i++)
218  {
219  vec_signal_samples.push_back(raw_tower->get_signal_samples(i));
220  }
221 
222  double peak = NAN;
223  double peak_sample = NAN;
224  double pedstal = NAN;
225 
226  switch (_fit_type)
227  {
228  case kPowerLawExp:
229  PROTOTYPE4_FEM::SampleFit_PowerLawExp(vec_signal_samples, peak,
230  peak_sample, pedstal, Verbosity());
231  break;
232 
233  case kPeakSample:
234  PROTOTYPE4_FEM::SampleFit_PeakSample(vec_signal_samples, peak,
235  peak_sample, pedstal, Verbosity());
236  break;
237 
238  case kPowerLawDoubleExp:
239  {
240  map<int, double> parameters_io;
241 
242  PROTOTYPE4_FEM::SampleFit_PowerLawDoubleExp(vec_signal_samples, peak,
243  peak_sample, pedstal,
244  parameters_io, Verbosity());
245  }
246  break;
247 
249  {
250  map<int, double> parameters_io(parameters_constraints);
251 
252  PROTOTYPE4_FEM::SampleFit_PowerLawDoubleExp(vec_signal_samples, peak,
253  peak_sample, pedstal,
254  parameters_io, Verbosity());
255  }
256  break;
257  default:
258  cout << __PRETTY_FUNCTION__ << " - FATAL error - unkown fit type "
259  << _fit_type << endl;
260  exit(3);
261  break;
262  }
263 
264  // store the result - raw_tower
265  if (std::isnan(raw_tower->get_energy()))
266  {
267  // Raw tower was never fit, store the current fit
268 
269  raw_tower->set_energy(peak);
270  raw_tower->set_time(peak_sample);
271 
272  if (Verbosity())
273  {
274  raw_tower->identify();
275  }
276  }
277 
278  // store the result - calib_tower
279  RawTower_Prototype4 *calib_tower = new RawTower_Prototype4(*raw_tower);
280  calib_tower->set_energy(peak * calibration_const);
281  calib_tower->set_time(peak_sample);
282 
283  for (int i = 0; i < RawTower_Prototype4::NSAMPLES; i++)
284  {
285  calib_tower->set_signal_samples(i, (vec_signal_samples[i] - pedstal) *
286  calibration_const);
287  }
288 
289  _calib_towers->AddTower(key, calib_tower);
290 
291  } // for (rtiter = begin_end.first; rtiter != begin_end.second; ++rtiter)
292 
293  if (Verbosity())
294  {
295  std::cout << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
296  << "input sum energy = " << _raw_towers->getTotalEdep()
297  << ", output sum digitalized value = "
298  << _calib_towers->getTotalEdep() << std::endl;
299  }
300 
302 }
303 
304 //_______________________________________
306 {
307  PHNodeIterator iter(topNode);
308 
309  PHCompositeNode *dstNode =
310  dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
311  if (!dstNode)
312  {
313  std::cerr << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
314  << "DST Node missing, doing nothing." << std::endl;
315  throw std::runtime_error(
316  "Failed to find DST node in RawTowerCalibration::CreateNodes");
317  }
318 
319  RawTowerNodeName = "TOWER_" + _raw_tower_node_prefix + "_" + detector;
320  _raw_towers =
321  findNode::getClass<RawTowerContainer>(dstNode, RawTowerNodeName.c_str());
322  if (!_raw_towers)
323  {
324  std::cerr << Name() << "::" << detector << "::" << __PRETTY_FUNCTION__
325  << " " << RawTowerNodeName << " Node missing, doing bail out!"
326  << std::endl;
327  throw std::runtime_error("Failed to find " + RawTowerNodeName +
328  " node in RawTowerCalibration::CreateNodes");
329  }
330 
331  // Create the tower nodes on the tree
332  PHNodeIterator dstiter(dstNode);
333  PHCompositeNode *DetNode = dynamic_cast<PHCompositeNode *>(
334  dstiter.findFirst("PHCompositeNode", detector));
335  if (!DetNode)
336  {
337  DetNode = new PHCompositeNode(detector);
338  dstNode->addNode(DetNode);
339  }
340 
341  // Be careful as a previous calibrator may have been registered for this
342  // detector
344  _calib_towers =
345  findNode::getClass<RawTowerContainer>(DetNode, CaliTowerNodeName.c_str());
346  if (!_calib_towers)
347  {
350  _calib_towers, CaliTowerNodeName.c_str(), "PHObject");
351  DetNode->addNode(towerNode);
352  }
353 
354  // update the parameters on the node tree
355  PHCompositeNode *parNode =
356  dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "RUN"));
357  assert(parNode);
358  const string paramnodename = string("Calibration_") + detector;
359 
360  // this step is moved to after detector construction
361  // save updated persistant copy on node tree
362  _calib_params.SaveToNodeTree(parNode, paramnodename);
363 }
364 
366 {
367  param.set_int_param("use_chan_calibration", 0);
368 
369  // additional scale for the calibration constant
370  // negative pulse -> positive with -1
371  param.set_double_param("calib_const_scale", 1);
372 }