CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CocoaDaqReaderRoot.cc
Go to the documentation of this file.
1 #include "../interface/CocoaDaqReaderRoot.h"
2 #include "TFile.h"
7 
9 
10 #include <iostream>
11 
12 #include "TClonesArray.h"
13 
14 
15 //----------------------------------------------------------------------
17 {
18  if ( ALIUtils::debug >= 3) 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;
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;
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 
52 //----------------------------------------------------------------------
54 {
55  theFile->Close();
56 }
57 
58 //----------------------------------------------------------------------
60 {
61  return ReadEvent( nextEvent );
62 }
63 
64 
65 //----------------------------------------------------------------------
67 {
68  std::vector<OpticalAlignMeasurementInfo> measList;
69 
70  int nb = 0; // dummy, number of bytes
71  // Loop over all events
72  nb = theTree->GetEntry(nev); // read in entire event
73 
74  if ( ALIUtils::debug >= 3) std::cout << "CocoaDaqReaderRoot reading event " << nev << " " << nb << std::endl;
75  if( nb == 0 ) return 0; //end of file reached??
76 
77  // Every n events, dump one to screen
78  int n = 1;
79  if(nev%n == 0 && ALIUtils::debug >= 3 ) theEvent->DumpIt();
80 
81  //if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent "<< nev <<std::endl;
82 
83  if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent npos2D "<< theEvent->GetNumPos2D() << " nCOPS " << theEvent->GetNumPosCOPS() << std::endl;
84 
85  for(int ii=0; ii<theEvent->GetNumPos2D(); ii++) {
87  if ( ALIUtils::debug >= 4) std::cout<<"2D sensor "<<ii<<" has ID = "<<pos2D->GetID()<< std::endl;
88  pos2D->DumpIt("2DSENSOR");
89  measList.push_back( GetMeasFromPosition2D( pos2D ) );
90  }
91  for(int ii=0; ii<theEvent->GetNumPosCOPS(); ii++) {
93  measList.push_back( GetMeasFromPositionCOPS( posCOPS ) );
94  if ( ALIUtils::debug >= 4) {
95  std::cout<<"COPS sensor "<<ii<<" has ID = "<<posCOPS->GetID()<< std::endl;
96  posCOPS->DumpIt("COPS");
97  }
98  }
99  for(int ii=0; ii<theEvent->GetNumTilt(); ii++) {
100  AliDaqTilt* tilt = (AliDaqTilt*) theEvent->GetArray_Tilt()->At(ii);
101  measList.push_back( GetMeasFromTilt( tilt ) );
102  if ( ALIUtils::debug >= 4) {
103  std::cout<<"TILT sensor "<<ii<<" has ID = "<<tilt->GetID()<< std::endl;
104  tilt->DumpIt("TILT");
105  }
106 
107  }
108  for(int ii=0; ii<theEvent->GetNumDist(); ii++) {
110  measList.push_back( GetMeasFromDist( dist ) );
111  if ( ALIUtils::debug >= 4) {
112  std::cout<<"DIST sensor "<<ii<<" has ID = "<<dist->GetID()<< std::endl;
113  dist->DumpIt("DIST");
114  }
115  }
116 
117  nextEvent = nev + 1;
118 
119  BuildMeasurementsFromOptAlign( measList );
120 
121  return 1;
122 
123 }
124 
125 //----------------------------------------------------------------------
127 {
129 
130  meas.type_ = "SENSOR2D";
131  meas.name_ = pos2D->GetID();
132  //- std::vector<std::string> measObjectNames_;
133  std::vector<bool> isSimu;
134  for( size_t jj = 0; jj < 2; jj++ ){
135  isSimu.push_back(false);
136  }
137  meas.isSimulatedValue_ = isSimu;
138  std::vector<OpticalAlignParam> paramList;
139  OpticalAlignParam oaParam1;
140  oaParam1.name_ = "H:";
141  oaParam1.value_ = pos2D->GetX()/100.;
142  oaParam1.error_ = pos2D->GetXerror()/100.;
143  paramList.push_back(oaParam1);
144 
145  OpticalAlignParam oaParam2;
146  oaParam2.name_ = "V:";
147  oaParam2.value_ = pos2D->GetY()/100.;
148  oaParam2.error_ = pos2D->GetYerror()/100.;
149  paramList.push_back(oaParam2);
150 
151  meas.values_ = paramList;
152 
153  return meas;
154 }
155 
156 
157 //----------------------------------------------------------------------
159 {
161 
162  meas.type_ = "COPS";
163  meas.name_ = posCOPS->GetID();
164  //- std::vector<std::string> measObjectNames_;
165  std::vector<bool> isSimu;
166  for( size_t jj = 0; jj < 4; jj++ ){
167  isSimu.push_back(false);
168  }
169  meas.isSimulatedValue_ = isSimu;
170 
171  std::vector<OpticalAlignParam> paramList;
172  OpticalAlignParam oaParam1;
173  oaParam1.name_ = "U:";
174  oaParam1.value_ = posCOPS->GetUp()/100.;
175  oaParam1.error_ = posCOPS->GetUpError()/100.;
176  paramList.push_back(oaParam1);
177 
178  OpticalAlignParam oaParam2;
179  oaParam2.name_ = "U:";
180  oaParam2.value_ = posCOPS->GetDown()/100.;
181  oaParam2.error_ = posCOPS->GetDownError()/100.;
182  paramList.push_back(oaParam2);
183 
184  OpticalAlignParam oaParam3;
185  oaParam3.name_ = "U:";
186  oaParam3.value_ = posCOPS->GetRight()/100.;
187  oaParam3.error_ = posCOPS->GetRightError()/100.;
188  paramList.push_back(oaParam3);
189 
190  OpticalAlignParam oaParam4;
191  oaParam4.name_ = "U:";
192  oaParam4.value_ = posCOPS->GetLeft()/100.;
193  oaParam4.error_ = posCOPS->GetLeftError()/100.;
194  paramList.push_back(oaParam4);
195 
196  meas.values_ = paramList;
197 
198  return meas;
199 
200 }
201 
202 //----------------------------------------------------------------------
204 {
206 
207  meas.type_ = "TILTMETER";
208  meas.name_ = tilt->GetID();
209  //- std::vector<std::string> measObjectNames_;
210  std::vector<bool> isSimu;
211  for( size_t jj = 0; jj < 2; jj++ ){
212  isSimu.push_back(false);
213  }
214  meas.isSimulatedValue_ = isSimu;
215  std::vector<OpticalAlignParam> paramList;
216  OpticalAlignParam oaParam;
217  oaParam.name_ = "T:";
218  oaParam.value_ = tilt->GetTilt();
219  oaParam.error_ = tilt->GetTiltError();
220  paramList.push_back(oaParam);
221 
222  meas.values_ = paramList;
223 
224  return meas;
225 
226 }
227 
228 
229 //----------------------------------------------------------------------
231 {
233 
234  meas.type_ = "DISTANCEMETER";
235  meas.name_ = dist->GetID();
236  //- std::vector<std::string> measObjectNames_;
237  std::vector<bool> isSimu;
238  for( size_t jj = 0; jj < 2; jj++ ){
239  isSimu.push_back(false);
240  }
241  meas.isSimulatedValue_ = isSimu;
242  std::vector<OpticalAlignParam> paramList;
243  OpticalAlignParam oaParam;
244  oaParam.name_ = "D:";
245  oaParam.value_ = dist->GetDistance()/100.;
246  oaParam.error_ = dist->GetDistanceError()/100.;
247  paramList.push_back(oaParam);
248 
249  meas.values_ = paramList;
250 
251  return meas;
252 }
253 
254 
255 //----------------------------------------------------------------------
256 void CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign( std::vector<OpticalAlignMeasurementInfo>& measList )
257 {
258  if ( ALIUtils::debug >= 3) std::cout << "@@@ CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign " << std::endl;
259 
260  //set date and time of current measurement
261  // if( wordlist[0] == "DATE:" ) {
262  // Measurement::setCurrentDate( wordlist );
263  // }
264 
265  //---------- loop measurements read from ROOT and check for corresponding measurement in Model
266  // ALIint nMeasModel = Model::MeasurementList().size();
267  ALIint nMeasRoot = measList.size();
268  if(ALIUtils::debug >= 4) {
269  std::cout << " Building " << nMeasRoot << " measurements from ROOT file " << std::endl;
270  }
271 
272  //--- Loop to Measurements in Model and check for corresponding measurement in ROOT
273  std::vector< Measurement* >::const_iterator vmcite;
274  for( vmcite = Model::MeasurementList().begin(); vmcite != Model::MeasurementList().end(); vmcite++ ) {
275  ALIint fcolon = (*vmcite)->name().find(':');
276  ALIstring oname = (*vmcite)->name();
277  oname = oname.substr(fcolon+1,oname.length());
278 
279  //---------- loop measurements read from ROOT
280  ALIint ii;
281  for(ii = 0; ii < nMeasRoot; ii++) {
282  OpticalAlignMeasurementInfo measInfo = measList[ii];
283  std::cout << " measurement name ROOT " << measInfo.name_ << " Model= " << (*vmcite)->name() << " short " << oname << std::endl;
284 
285  if( oname == measInfo.name_ ) {
286  //-------- Measurement found, fill data
287  //---- Check that type is the same
288  if( (*vmcite)->type() != measInfo.type_ ) {
289  std::cerr << "!!! Measurement from ROOT file: type in file is "
290  <<measInfo.type_ << " and should be " << (*vmcite)->type() << std::endl;
291  exit(1);
292  }
293 
294  std::cout << " NOBJECTS IN MEAS " << (*vmcite)->OptOList().size() << " NMEAS " << Model::MeasurementList().size() << std::endl;
295 
296  std::vector<OpticalAlignParam> measValues = measInfo.values_;
297 
298  for( size_t jj= 0; jj < measValues.size(); jj++ ){
299  (*vmcite)->fillData( jj, &(measValues[jj]) );
300  }
301 
302  std::cout << " NOBJECTS IN MEAS after " << (*vmcite)->OptOList().size() << " NMEAS " << Model::MeasurementList().size() << std::endl;
303 
304  break;
305  }
306  }
307  if (ii==nMeasRoot) {
308  std::cerr << "!!! Reading measurement from file: measurement not found! Type in list is " << oname << std::endl;
309  exit(1);
310  }
311  }
312 
313 }
314 
Float_t GetX() const
Float_t GetTilt() const
void DumpIt(TString Name)
TClonesArray * GetArray_Tilt() const
OpticalAlignMeasurementInfo GetMeasFromTilt(AliDaqTilt *tilt)
OpticalAlignMeasurementInfo GetMeasFromDist(AliDaqDistance *dist)
int ALIint
Definition: CocoaGlobals.h:15
TClonesArray * GetArray_PositionCOPS() const
static ALIint debug
Definition: ALIUtils.h:35
Float_t GetUpError() const
int ii
Definition: cuy.py:588
TString GetID()
TClonesArray * GetArray_Position2D() const
Float_t GetYerror() const
CocoaDaqReaderRoot(const std::string &m_inFileName)
int GetNumTilt() const
virtual void BuildMeasurementsFromOptAlign(std::vector< OpticalAlignMeasurementInfo > &measList)
CocoaDaqRootEvent * theEvent
virtual bool ReadNextEvent()
Float_t GetDistanceError() const
Float_t GetLeft() const
Float_t GetRight() const
int GetNumPos2D() const
Float_t GetTiltError() const
void DumpIt(TString Name)
void DumpIt(TString Name)
OpticalAlignMeasurementInfo GetMeasFromPositionCOPS(AliDaqPositionCOPS *posCOPS)
TClonesArray * GetArray_Dist() const
static void SetDaqReader(CocoaDaqReader *reader)
virtual bool ReadEvent(int nev)
Float_t GetXerror() const
Float_t GetY() const
Float_t GetLeftError() const
int GetNumDist() const
Float_t GetDistance() const
#define begin
Definition: vmac.h:31
Float_t GetDownError() const
std::string ALIstring
Definition: CocoaGlobals.h:9
std::vector< OpticalAlignParam > values_
void DumpIt(TString Name)
tuple cout
Definition: gather_cfg.py:121
OpticalAlignMeasurementInfo GetMeasFromPosition2D(AliDaqPosition2D *pos2D)
Float_t GetDown() const
Float_t GetUp() const
static std::vector< Measurement * > & MeasurementList()
Definition: Model.h:79
Float_t GetRightError() const
int GetNumPosCOPS() const