CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CocoaDaqReaderRoot.cc
Go to the documentation of this file.
2 #include "TFile.h"
7 
9 
10 #include <iostream>
11 #include <string>
12 
13 #include "TClonesArray.h"
14 
15 //----------------------------------------------------------------------
17  if (ALIUtils::debug >= 3)
18  std::cout << " CocoaDaqReaderRoot opening file: " << m_inFileName << std::endl;
19  // Open root file
20  theFile = new TFile(m_inFileName.c_str());
21  if (!theTree) {
22  std::cerr << " CocoaDaqReaderRoot TTree file not found " << m_inFileName << std::endl;
23  throw std::exception();
24  }
25 
26  // Read TTree named "CocoaDaq" in memory. !! SHOULD BE CALLED Alignment_Cocoa
27  theTree = (TTree*)theFile->Get("CocoaDaq");
28  // theTree = (TTree*)theFile->Get("Alignment_Link_Cocoa");
29 
30  if (!theTree) {
31  std::cerr << " CocoaDaqReaderRoot TTree in file " << m_inFileName << " should be called 'CocoaDaq' " << std::endl;
32  throw std::exception();
33  }
34  TBranch* branch = theTree->GetBranch("Alignment_Cocoa");
35 
36  nev = branch->GetEntries(); // number of entries in Tree
37  //if ( ALIUtils::debug >= 2) std::cout << "CocoaDaqReaderRoot::CocoaDaqReaderRoot: number of entries in Tree " << nev << std::endl;
38 
39  nextEvent = 0;
40 
41  // Event object must be created before setting the branch address
43 
44  // link pointer to Tree branch
45  theTree->SetBranchAddress("Alignment_Cocoa", &theEvent); // !! SHOULD BE CALLED Alignment_Cocoa
46  // theTree->SetBranchAddress("Alignment_Link", &theEvent); // !! SHOULD BE CALLED Alignment_Cocoa
47 
49 }
50 
51 //----------------------------------------------------------------------
53 
54 //----------------------------------------------------------------------
56 
57 //----------------------------------------------------------------------
59  std::vector<OpticalAlignMeasurementInfo> measList;
60 
61  int nb = 0; // dummy, number of bytes
62  // Loop over all events
63  nb = theTree->GetEntry(nev); // read in entire event
64 
65  if (ALIUtils::debug >= 3)
66  std::cout << "CocoaDaqReaderRoot reading event " << nev << " " << nb << std::endl;
67  if (nb == 0)
68  return false; //end of file reached??
69 
70  // Every n events, dump one to screen
71  int n = 1;
72  if (nev % n == 0 && ALIUtils::debug >= 3)
73  theEvent->DumpIt();
74 
75  //if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent "<< nev <<std::endl;
76 
77  if (ALIUtils::debug >= 3)
78  std::cout << " CocoaDaqReaderRoot::ReadEvent npos2D " << theEvent->GetNumPos2D() << " nCOPS "
79  << theEvent->GetNumPosCOPS() << std::endl;
80 
81  for (int ii = 0; ii < theEvent->GetNumPos2D(); ii++) {
83  if (ALIUtils::debug >= 4)
84  std::cout << "2D sensor " << ii << " has ID = " << pos2D->GetID() << std::endl;
85  pos2D->DumpIt("2DSENSOR");
86  measList.push_back(GetMeasFromPosition2D(pos2D));
87  }
88  for (int ii = 0; ii < theEvent->GetNumPosCOPS(); ii++) {
90  measList.push_back(GetMeasFromPositionCOPS(posCOPS));
91  if (ALIUtils::debug >= 4) {
92  std::cout << "COPS sensor " << ii << " has ID = " << posCOPS->GetID() << std::endl;
93  posCOPS->DumpIt("COPS");
94  }
95  }
96  for (int ii = 0; ii < theEvent->GetNumTilt(); ii++) {
97  AliDaqTilt* tilt = (AliDaqTilt*)theEvent->GetArray_Tilt()->At(ii);
98  measList.push_back(GetMeasFromTilt(tilt));
99  if (ALIUtils::debug >= 4) {
100  std::cout << "TILT sensor " << ii << " has ID = " << tilt->GetID() << std::endl;
101  tilt->DumpIt("TILT");
102  }
103  }
104  for (int ii = 0; ii < theEvent->GetNumDist(); ii++) {
106  measList.push_back(GetMeasFromDist(dist));
107  if (ALIUtils::debug >= 4) {
108  std::cout << "DIST sensor " << ii << " has ID = " << dist->GetID() << std::endl;
109  dist->DumpIt("DIST");
110  }
111  }
112 
113  nextEvent = nev + 1;
114 
116 
117  return true;
118 }
119 
120 //----------------------------------------------------------------------
123 
124  meas.type_ = "SENSOR2D";
125  meas.name_ = std::string(pos2D->GetID().Data());
126  //- std::vector<std::string> measObjectNames_;
127  std::vector<bool> isSimu;
128  for (size_t jj = 0; jj < 2; jj++) {
129  isSimu.push_back(false);
130  }
131  meas.isSimulatedValue_ = isSimu;
132  std::vector<OpticalAlignParam> paramList;
133  OpticalAlignParam oaParam1;
134  oaParam1.name_ = "H:";
135  oaParam1.value_ = pos2D->GetX() / 100.;
136  oaParam1.error_ = pos2D->GetXerror() / 100.;
137  paramList.push_back(oaParam1);
138 
139  OpticalAlignParam oaParam2;
140  oaParam2.name_ = "V:";
141  oaParam2.value_ = pos2D->GetY() / 100.;
142  oaParam2.error_ = pos2D->GetYerror() / 100.;
143  paramList.push_back(oaParam2);
144 
145  meas.values_ = paramList;
146 
147  return meas;
148 }
149 
150 //----------------------------------------------------------------------
153 
154  meas.type_ = "COPS";
155  meas.name_ = std::string(posCOPS->GetID().Data());
156  //- std::vector<std::string> measObjectNames_;
157  std::vector<bool> isSimu;
158  for (size_t jj = 0; jj < 4; jj++) {
159  isSimu.push_back(false);
160  }
161  meas.isSimulatedValue_ = isSimu;
162 
163  std::vector<OpticalAlignParam> paramList;
164  OpticalAlignParam oaParam1;
165  oaParam1.name_ = "U:";
166  oaParam1.value_ = posCOPS->GetUp() / 100.;
167  oaParam1.error_ = posCOPS->GetUpError() / 100.;
168  paramList.push_back(oaParam1);
169 
170  OpticalAlignParam oaParam2;
171  oaParam2.name_ = "U:";
172  oaParam2.value_ = posCOPS->GetDown() / 100.;
173  oaParam2.error_ = posCOPS->GetDownError() / 100.;
174  paramList.push_back(oaParam2);
175 
176  OpticalAlignParam oaParam3;
177  oaParam3.name_ = "U:";
178  oaParam3.value_ = posCOPS->GetRight() / 100.;
179  oaParam3.error_ = posCOPS->GetRightError() / 100.;
180  paramList.push_back(oaParam3);
181 
182  OpticalAlignParam oaParam4;
183  oaParam4.name_ = "U:";
184  oaParam4.value_ = posCOPS->GetLeft() / 100.;
185  oaParam4.error_ = posCOPS->GetLeftError() / 100.;
186  paramList.push_back(oaParam4);
187 
188  meas.values_ = paramList;
189 
190  return meas;
191 }
192 
193 //----------------------------------------------------------------------
196 
197  meas.type_ = "TILTMETER";
198  meas.name_ = std::string(tilt->GetID().Data());
199  //- std::vector<std::string> measObjectNames_;
200  std::vector<bool> isSimu;
201  for (size_t jj = 0; jj < 2; jj++) {
202  isSimu.push_back(false);
203  }
204  meas.isSimulatedValue_ = isSimu;
205  std::vector<OpticalAlignParam> paramList;
206  OpticalAlignParam oaParam;
207  oaParam.name_ = "T:";
208  oaParam.value_ = tilt->GetTilt();
209  oaParam.error_ = tilt->GetTiltError();
210  paramList.push_back(oaParam);
211 
212  meas.values_ = paramList;
213 
214  return meas;
215 }
216 
217 //----------------------------------------------------------------------
220 
221  meas.type_ = "DISTANCEMETER";
222  meas.name_ = std::string(dist->GetID().Data());
223  //- std::vector<std::string> measObjectNames_;
224  std::vector<bool> isSimu;
225  for (size_t jj = 0; jj < 2; jj++) {
226  isSimu.push_back(false);
227  }
228  meas.isSimulatedValue_ = isSimu;
229  std::vector<OpticalAlignParam> paramList;
230  OpticalAlignParam oaParam;
231  oaParam.name_ = "D:";
232  oaParam.value_ = dist->GetDistance() / 100.;
233  oaParam.error_ = dist->GetDistanceError() / 100.;
234  paramList.push_back(oaParam);
235 
236  meas.values_ = paramList;
237 
238  return meas;
239 }
240 
241 //----------------------------------------------------------------------
242 void CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign(std::vector<OpticalAlignMeasurementInfo>& measList) {
243  if (ALIUtils::debug >= 3)
244  std::cout << "@@@ CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign " << std::endl;
245 
246  //set date and time of current measurement
247  // if( wordlist[0] == "DATE:" ) {
248  // Measurement::setCurrentDate( wordlist );
249  // }
250 
251  //---------- loop measurements read from ROOT and check for corresponding measurement in Model
252  // ALIint nMeasModel = Model::MeasurementList().size();
253  ALIint nMeasRoot = measList.size();
254  if (ALIUtils::debug >= 4) {
255  std::cout << " Building " << nMeasRoot << " measurements from ROOT file " << std::endl;
256  }
257 
258  //--- Loop to Measurements in Model and check for corresponding measurement in ROOT
259  std::vector<Measurement*>::const_iterator vmcite;
260  for (vmcite = Model::MeasurementList().begin(); vmcite != Model::MeasurementList().end(); ++vmcite) {
261  ALIint fcolon = (*vmcite)->name().find(':');
262  ALIstring oname = (*vmcite)->name();
263  oname = oname.substr(fcolon + 1, oname.length());
264 
265  //---------- loop measurements read from ROOT
266  ALIint ii;
267  for (ii = 0; ii < nMeasRoot; ii++) {
268  OpticalAlignMeasurementInfo measInfo = measList[ii];
269  std::cout << " measurement name ROOT " << measInfo.name_ << " Model= " << (*vmcite)->name() << " short " << oname
270  << std::endl;
271 
272  if (oname == measInfo.name_) {
273  //-------- Measurement found, fill data
274  //---- Check that type is the same
275  if ((*vmcite)->type() != measInfo.type_) {
276  std::cerr << "!!! Measurement from ROOT file: type in file is " << measInfo.type_ << " and should be "
277  << (*vmcite)->type() << std::endl;
278  exit(1);
279  }
280 
281  std::cout << " NOBJECTS IN MEAS " << (*vmcite)->OptOList().size() << " NMEAS "
282  << Model::MeasurementList().size() << std::endl;
283 
284  std::vector<OpticalAlignParam> measValues = measInfo.values_;
285 
286  for (size_t jj = 0; jj < measValues.size(); jj++) {
287  (*vmcite)->fillData(jj, &(measValues[jj]));
288  }
289 
290  std::cout << " NOBJECTS IN MEAS after " << (*vmcite)->OptOList().size() << " NMEAS "
291  << Model::MeasurementList().size() << std::endl;
292 
293  break;
294  }
295  }
296  if (ii == nMeasRoot) {
297  std::cerr << "!!! Reading measurement from file: measurement not found! Type in list is " << oname << std::endl;
298  exit(1);
299  }
300  }
301 }
Float_t GetX() const
void DumpIt(const TString &Name)
Float_t GetTilt() const
TClonesArray * GetArray_Tilt() const
OpticalAlignMeasurementInfo GetMeasFromTilt(AliDaqTilt *tilt)
bool ReadEvent(int nev) override
OpticalAlignMeasurementInfo GetMeasFromDist(AliDaqDistance *dist)
int ALIint
Definition: CocoaGlobals.h:15
TClonesArray * GetArray_PositionCOPS() const
static ALIint debug
Definition: ALIUtils.h:34
~CocoaDaqReaderRoot() override
Float_t GetUpError() const
int ii
Definition: cuy.py:589
TString GetID()
TClonesArray * GetArray_Position2D() const
Float_t GetYerror() const
CocoaDaqReaderRoot(const std::string &m_inFileName)
int GetNumTilt() const
CocoaDaqRootEvent * theEvent
Float_t GetDistanceError() const
bool ReadNextEvent() override
Float_t GetLeft() const
Float_t GetRight() const
int GetNumPos2D() const
Float_t GetTiltError() const
void DumpIt(const TString &Name)
OpticalAlignMeasurementInfo GetMeasFromPositionCOPS(AliDaqPositionCOPS *posCOPS)
TClonesArray * GetArray_Dist() const
static void SetDaqReader(CocoaDaqReader *reader)
Float_t GetXerror() const
Float_t GetY() const
Float_t GetLeftError() const
string oname
Definition: heppy_report.py:57
int GetNumDist() const
Float_t GetDistance() const
void BuildMeasurementsFromOptAlign(std::vector< OpticalAlignMeasurementInfo > &measList) override
Float_t GetDownError() const
std::string ALIstring
Definition: CocoaGlobals.h:9
std::vector< OpticalAlignParam > values_
tuple cout
Definition: gather_cfg.py:144
void DumpIt(const TString &Name)
static std::vector< Measurement * > & MeasurementList()
Definition: Model.h:88
OpticalAlignMeasurementInfo GetMeasFromPosition2D(AliDaqPosition2D *pos2D)
Float_t GetDown() const
Float_t GetUp() const
Float_t GetRightError() const
void DumpIt(const TString &Name)
int GetNumPosCOPS() const