Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DDIRCLutReader.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DDIRCLutReader.cc
4 //
5 
6 #include <unistd.h>
7 #include <cassert>
8 #include <math.h>
9 using namespace std;
10 
11 #include "DDIRCLutReader.h"
12 #include "DANA/DApplication.h"
13 
14 //---------------------------------
15 // DDIRCLutReader (Constructor)
16 //---------------------------------
17 DDIRCLutReader::DDIRCLutReader(JApplication *japp, unsigned int run_number)
18 {
19  /////////////////////////////////////////
20  // retrieve from LUT from file or CCDB //
21  /////////////////////////////////////////
22  const int luts = 48;
23 
24  string lut_file;
25  gPARMS->SetDefaultParameter("DIRC_LUT", lut_file, "DIRC LUT root file (will eventually be moved to resource)");
26 
27  // follow similar procedure as other resources (DMagneticFieldMapFineMesh)
28  map<string,string> lut_map_name;
29  jcalib = japp->GetJCalibration(run_number);
30  if(jcalib->GetCalib("/DIRC/LUT/lut_map", lut_map_name))
31  jout << "Can't find requested /DIRC/LUT/lut_map in CCDB for this run!" << endl;
32  else if(lut_map_name.find("map_name") != lut_map_name.end() && lut_map_name["map_name"] != "None") {
33  jresman = japp->GetJResourceManager(run_number);
34  lut_file = jresman->GetResource(lut_map_name["map_name"]);
35  }
36 
37  // only create reader if have tree from CCDB or command-line parameter is set
38  if(!lut_file.empty()) {
39  jout<<"Reading DIRC LUT TTree from "<<lut_file<<" ..."<<endl;
40 
41  auto saveDir = gDirectory;
42  TFile *fLut = new TFile(lut_file.c_str());
43  if( !fLut->IsOpen() ){
44  jerr << "Unable to open " << lut_file << "!!" << endl;
45  _exit(-1);
46  }
47  TTree *tLut=(TTree*) fLut->Get("lut_dirc_flat");
48  if( tLut == NULL ){
49  jerr << "Unable find TTree lut_dirc_flat in " << lut_file << "!!" << endl;
50  _exit(-1);
51  }
52 
53  vector<Float_t> *LutPixelAngleX[luts];
54  vector<Float_t> *LutPixelAngleY[luts];
55  vector<Float_t> *LutPixelAngleZ[luts];
56  vector<Float_t> *LutPixelTime[luts];
57  vector<Long64_t> *LutPixelPath[luts];
58 
59  // clear arrays to fill from TTree
60  for(int l=0; l<luts; l++){
61  LutPixelAngleX[l] = 0;
62  LutPixelAngleY[l] = 0;
63  LutPixelAngleZ[l] = 0;
64  LutPixelTime[l] = 0;
65  LutPixelPath[l] = 0;
66  }
67 
68  for(int l=0; l<luts; l++){
69  tLut->SetBranchAddress(Form("LUT_AngleX_%d",l),&LutPixelAngleX[l]);
70  tLut->SetBranchAddress(Form("LUT_AngleY_%d",l),&LutPixelAngleY[l]);
71  tLut->SetBranchAddress(Form("LUT_AngleZ_%d",l),&LutPixelAngleZ[l]);
72  tLut->SetBranchAddress(Form("LUT_Time_%d",l),&LutPixelTime[l]);
73  tLut->SetBranchAddress(Form("LUT_Path_%d",l),&LutPixelPath[l]);
74  }
75 
76  // fill nodes with LUT info for each bar/pixel combination
77  for(int i=0; i<tLut->GetEntries(); i++) { // get pixels from TTree
78  tLut->GetEntry(i);
79 
80  for(int l=0; l<luts; l++){ // loop over bars
81  for(uint j=0; j<LutPixelAngleX[l]->size(); j++) { // loop over possible paths
82  TVector3 angle(LutPixelAngleX[l]->at(j), LutPixelAngleY[l]->at(j), LutPixelAngleZ[l]->at(j));
83  lutNodeAngle[l][i].push_back(angle);
84  lutNodeTime[l][i].push_back(LutPixelTime[l]->at(j));
85  lutNodePath[l][i].push_back(LutPixelPath[l]->at(j));
86  }
87  }
88  }
89 
90  // close LUT file
91  fLut->Close();
92  saveDir->cd();
93  }
94 }
95 
97 
98 }
99 
100 uint DDIRCLutReader::GetLutPixelAngleSize(int bar, int pixel) const
101 {
102 
103  return lutNodeAngle[bar][pixel].size();
104 }
105 
106 uint DDIRCLutReader::GetLutPixelTimeSize(int bar, int pixel) const
107 {
108  return lutNodeTime[bar][pixel].size();
109 }
110 
111 uint DDIRCLutReader::GetLutPixelPathSize(int bar, int pixel) const
112 {
113  return lutNodePath[bar][pixel].size();
114 }
115 
116 TVector3 DDIRCLutReader::GetLutPixelAngle(int bar, int pixel, int entry) const
117 {
118  return lutNodeAngle[bar][pixel].at(entry);
119 }
120 
121 Float_t DDIRCLutReader::GetLutPixelTime(int bar, int pixel, int entry) const
122 {
123  return lutNodeTime[bar][pixel].at(entry);
124 }
125 
126 Long64_t DDIRCLutReader::GetLutPixelPath(int bar, int pixel, int entry) const
127 {
128  return lutNodePath[bar][pixel].at(entry);
129 }
TVector3 GetLutPixelAngle(int bar, int pixel, int entry) const
uint GetLutPixelAngleSize(int bar, int pixel) const
Long64_t GetLutPixelPath(int bar, int pixel, int entry) const
Float_t GetLutPixelTime(int bar, int pixel, int entry) const
JApplication * japp
uint GetLutPixelPathSize(int bar, int pixel) const
virtual ~DDIRCLutReader()
DDIRCLutReader(JApplication *japp, unsigned int run_number)
uint GetLutPixelTimeSize(int bar, int pixel) const