Belle II Software  release-06-01-15
eclComputePulseTemplates_Step0.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include<iostream>
10 #include<fstream>
11 #include<stdlib.h>
12 #include"TTree.h"
13 #include"TString.h"
14 #include"TFile.h"
15 #include <assert.h>
16 
17 using namespace std;
18 
19 /*
20 Note photon shape calibrations are needed to compute hadron and diode templates.
21 Photon calibrations should be placed in file named params_gamma_shape.dat in same directiory as this script.
22 Format of params_gamma_shape.dat:
23 
24 cellID PhotonScale 0 p0 p1 p2 p3 p4 p5 p6 p7 p8 p9
25 
26 where p0-p9 are 10 ShaperDSP shaper parameters for basf2 C++ implementation of ShaperDSP function.
27 
28 eclComputePulseTemplates_Step0.cc - places photon parameters in ntuple
29 eclComputePulseTemplates_Step1.cc - produce photon array input for eclComputePulseTemplates_Step2.py
30 eclComputePulseTemplates_Step2.py - produce hadron and diode array input for eclComputePulseTemplates_Step3.py
31 eclComputePulseTemplates_Step3.cc - computed hadron and diode shape parameters and organizes into temp ntuple
32 eclComputePulseTemplates_Step4.cc - organizes photon, hadron and diode templates into final ntuple to produce dbobject.
33 
34 To run template calculation:
35 
36 First "OutputDirectory" variable in eclComputePulseTemplates_Stepx.cc files should be modified to a temp directory before running scripts.
37 recompile with "scons ecl"
38 
39 1) execute eclComputePulseTemplates_Step0 0
40 2) for i=0 i<873 i++
41 3) execute eclComputePulseTemplates_Step1 i*10 (i+1)*10
42 4) basf2 eclComputePulseTemplates_Step2.py i*10 (i+1)*10
43 5) execute eclComputePulseTemplates_Step3 i*10 (i+1)*10
44 6) execute eclComputePulseTemplates_Step0 1
45 7) execute eclWriteWaveformParametersLocalDB
46 Final step will create local db object.
47 */
48 
50 struct crystalInfo {
51  vector<double> PhotonWaveformPars;
52  vector<double> HadronWaveformPars;
53  vector<double> DiodeWaveformPars;
54  double MaxResHadron;
55  double MaxValHadron;
56  double MaxResDiode;
57  double MaxValDiode;
58 };
59 
60 int main(int argc, char* argv[])
61 {
62  //
63  TString OutputDirectory = "";
64  if (OutputDirectory == "") {
65  std::cout << "Error set ouput directory" << std::endl;
66  return -1;
67  }
68  //
69  assert(argc == 2);
70  int Flag = atoi(argv[1]);
71  TString InputDirectory = OutputDirectory;
72  //
73  TString ParameterTreeOutputName = OutputDirectory + "PhotonWaveformParameters.root";
74  if (Flag == 1) ParameterTreeOutputName = "DigitWaveformParameters.root";
75  TFile* ParameterTreeOutputFile = new TFile(ParameterTreeOutputName, "RECREATE");
76  //
77  TTree* ParameterTree = new TTree("ParTree", "");
78  double PhotonWaveformPar[11];
79  double HadronWaveformPar[11];
80  double DiodeWaveformPar[11];
81  for (int k = 0; k < 11; k++) {
82  PhotonWaveformPar[k] = 0;
83  HadronWaveformPar[k] = 0;
84  DiodeWaveformPar[k] = 0;
85  }
86  double mHadronRes = 0;
87  double mDiodeRes = 0;
88  double mHadronMax = 0;
89  double mDiodeMax = 0;
90  ParameterTree->Branch("PhotonPar", &PhotonWaveformPar, "PhotonWaveformPar[11]/D");
91  ParameterTree->Branch("HadronPar", &HadronWaveformPar, "HadronWaveformPar[11]/D");
92  ParameterTree->Branch("DiodePar", &DiodeWaveformPar, "DiodeWaveformPar[11]/D");
93  ParameterTree->Branch("mHadronRes", &mHadronRes, "mHadronRes/D");
94  ParameterTree->Branch("mDiodeRes", &mDiodeRes, "mDiodeRes/D");
95  ParameterTree->Branch("mHadronMax", &mHadronMax, "mHadronMax/D");
96  ParameterTree->Branch("mDiodeMax", &mDiodeMax, "mDiodeMax/D");
97  //
98  std::vector<crystalInfo> cellIDcheck(8736);
99  ifstream PhotonFile("/home/belle2/longos/WaveformFitting/ecl/tools/params_gamma_shape.dat");
100  if (PhotonFile.is_open()) {
101  std::vector<double> templine(12);
102  for (int k = 0; k < 8736; k++) {
103  for (unsigned int j = 0; j < templine.size(); j++) PhotonFile >> templine[j];
104  cellIDcheck[int(templine[0]) - 1].PhotonWaveformPars.resize(11);
105  cellIDcheck[int(templine[0]) - 1].PhotonWaveformPars[0] = templine[1];
106  for (int j = 0; j < 10; j++) cellIDcheck[int(templine[0]) - 1].PhotonWaveformPars[j + 1] = templine[j + 2];
107  std::cout << int(templine[0]) << " " << templine[1] << endl;
108  }
109  PhotonFile.close();
110  } else {
111  std::cout << "ERROR cannont open photon param file." << endl;
112  return -1;
113  }
114  if (Flag == 0) {
115  for (unsigned int f = 0; f < cellIDcheck.size(); f++) {
116  for (int k = 0; k < 11; k++) PhotonWaveformPar[k] = cellIDcheck[f].PhotonWaveformPars[k];
117  ParameterTree->Fill();
118  }
119  } else if (Flag == 1) {
120  int batch = 10;
121  for (int k = 0; k < 874; k++) {
122  int low = k * batch;
123  int high = (k + 1) * batch;
124  //std::cout<<k<<std::endl;
125  TFile* TempFile = new TFile(InputDirectory + Form("HadronPars_Low%d_High%d.root", low, high), "READ");
126  TTree* TempTree = (TTree*) TempFile->Get("HadronWaveformInfo");
127  double tHadronShapePars_A[11];
128  double tDiodeShapePars_A[11];
129  double tMaxResidualHadron_A;
130  double tMaxResidualDiode_A;
131  double tMaxValDiode_A;
132  double tMaxValHadron_A;
133  TempTree->SetBranchAddress("TempHadronPar11_A", &tHadronShapePars_A);
134  TempTree->SetBranchAddress("TempDiodePar11_A", &tDiodeShapePars_A);
135  TempTree->SetBranchAddress("MaxResHadron_A", &tMaxResidualHadron_A);
136  TempTree->SetBranchAddress("MaxResDiode_A", &tMaxResidualDiode_A);
137  TempTree->SetBranchAddress("MaxValDiode_A", &tMaxValDiode_A);
138  TempTree->SetBranchAddress("MaxValHadron_A", &tMaxValHadron_A);
139  //
140  for (int j = 0; j < batch; j++) {
141  int tCellID = (k * batch) + j;
142  if (tCellID >= 8736)continue;
143  TempTree->GetEntry(j);
144  cellIDcheck[tCellID].HadronWaveformPars.resize(11);
145  cellIDcheck[tCellID].DiodeWaveformPars.resize(11);
146  for (int p = 0; p < 11; p++) {
147  if (tHadronShapePars_A[p] > 100 || tHadronShapePars_A[p] < -100) {
148  std::cout << "Warning large parameter for: " << tCellID << " " << tHadronShapePars_A[p] << std::endl;
149  for (int h = 0; h < 11; h++)std::cout << tHadronShapePars_A[h] << " , ";
150  std::cout << std::endl;
151  }
152  cellIDcheck[tCellID].HadronWaveformPars[p] = tHadronShapePars_A[p];
153  }
154  //std::cout<<tCellID<<std::endl;
155  //for(int h=0;h<11;h++ )std::cout<<tHadronShapePars_A[h]<<" , ";
156  //std::cout<<std::endl;
157  cellIDcheck[tCellID].MaxResHadron = tMaxResidualHadron_A;
158  cellIDcheck[tCellID].MaxValHadron = tMaxValHadron_A;
159  for (int p = 0; p < 11; p++) cellIDcheck[tCellID].DiodeWaveformPars[p] = tDiodeShapePars_A[p];
160  cellIDcheck[tCellID].MaxResDiode = tMaxResidualDiode_A;
161  cellIDcheck[tCellID].MaxValDiode = tMaxValDiode_A;
162  }
163  TempFile->Close();
164  }
165  //
166  for (unsigned int f = 0; f < 8736; f++) {
167  for (int k = 0; k < 11; k++) {
168  PhotonWaveformPar[k] = cellIDcheck[f].PhotonWaveformPars[k];
169  HadronWaveformPar[k] = cellIDcheck[f].HadronWaveformPars[k];
170  DiodeWaveformPar[k] = cellIDcheck[f].DiodeWaveformPars[k];
171  }
172  mHadronRes = cellIDcheck[f].MaxResHadron;
173  mDiodeRes = cellIDcheck[f].MaxResDiode;
174  mHadronMax = cellIDcheck[f].MaxValHadron;
175  mDiodeMax = cellIDcheck[f].MaxValDiode;
176  ParameterTree->Fill();
177  }
178  }
179  //
180  ParameterTreeOutputFile->cd();
181  ParameterTree->Write();
182  ParameterTreeOutputFile->Write();
183  ParameterTreeOutputFile->Close();
184  //
185  //
186  return 0;
187 }
a struct to hold relavent Crystal information
vector< double > PhotonWaveformPars
photon waveform parameters
vector< double > HadronWaveformPars
hadron waveform parameters
vector< double > DiodeWaveformPars
diode waveform parameters
int main(int argc, char **argv)
Run all tests.
Definition: test_main.cc:75