CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/CocoaModel/src/CocoaDaqReaderRoot.cc

Go to the documentation of this file.
00001 #include "../interface/CocoaDaqReaderRoot.h"
00002 #include "TFile.h" 
00003 #include "Alignment/CocoaDaq/interface/CocoaDaqRootEvent.h"
00004 #include "Alignment/CocoaModel/interface/Measurement.h"
00005 #include "Alignment/CocoaModel/interface/Model.h"
00006 #include "Alignment/CocoaUtilities/interface/ALIUtils.h"
00007 
00008 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurements.h"
00009 
00010 #include <iostream>
00011 
00012 #include "TClonesArray.h"
00013 
00014 
00015 //----------------------------------------------------------------------
00016 CocoaDaqReaderRoot::CocoaDaqReaderRoot(const std::string& m_inFileName )
00017 {
00018   if ( ALIUtils::debug >= 3) std::cout << " CocoaDaqReaderRoot opening file: " << m_inFileName << std::endl;
00019   // Open root file
00020   theFile = new TFile(m_inFileName.c_str()); 
00021   if( !theTree ) {
00022     std::cerr << " CocoaDaqReaderRoot TTree file not found " << m_inFileName << std::endl;
00023     std::exception();
00024   }
00025   
00026   // Read TTree named "CocoaDaq" in memory.  !! SHOULD BE CALLED Alignment_Cocoa
00027    theTree = (TTree*)theFile->Get("CocoaDaq");
00028   //  theTree = (TTree*)theFile->Get("Alignment_Link_Cocoa");
00029   
00030   if( !theTree ) {
00031     std::cerr << " CocoaDaqReaderRoot TTree in file " << m_inFileName << " should be called 'CocoaDaq' " << std::endl;
00032     std::exception();
00033   }
00034   TBranch *branch = theTree->GetBranch("Alignment_Cocoa");
00035 
00036   nev = branch->GetEntries(); // number of entries in Tree
00037   //if ( ALIUtils::debug >= 2) std::cout << "CocoaDaqReaderRoot::CocoaDaqReaderRoot:  number of entries in Tree " << nev << std::endl;
00038  
00039   nextEvent = 0;
00040 
00041   // Event object must be created before setting the branch address
00042   theEvent = new CocoaDaqRootEvent();
00043 
00044   // link pointer to Tree branch
00045    theTree->SetBranchAddress("Alignment_Cocoa", &theEvent);  //  !! SHOULD BE CALLED Alignment_Cocoa
00046   // theTree->SetBranchAddress("Alignment_Link", &theEvent);  //  !! SHOULD BE CALLED Alignment_Cocoa
00047 
00048   CocoaDaqReader::SetDaqReader( this );
00049 
00050 }
00051 
00052 //----------------------------------------------------------------------
00053 CocoaDaqReaderRoot::~CocoaDaqReaderRoot()
00054 {
00055   theFile->Close();
00056 }
00057 
00058 //----------------------------------------------------------------------
00059 bool CocoaDaqReaderRoot::ReadNextEvent()
00060 {
00061   return ReadEvent( nextEvent );
00062 }
00063 
00064 
00065 //----------------------------------------------------------------------
00066 bool CocoaDaqReaderRoot::ReadEvent( int nev )
00067 {
00068   std::vector<OpticalAlignMeasurementInfo> measList;
00069 
00070   int nb  = 0;   // dummy, number of bytes
00071   // Loop over all events
00072   nb = theTree->GetEntry(nev);  // read in entire event
00073  
00074   if ( ALIUtils::debug >= 3) std::cout << "CocoaDaqReaderRoot reading event " << nev << " " << nb << std::endl;
00075   if( nb == 0 ) return 0; //end of file reached??
00076 
00077   // Every n events, dump one to screen
00078   int n = 1;
00079   if(nev%n == 0 &&  ALIUtils::debug >= 3 ) theEvent->DumpIt();
00080   
00081   //if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent "<< nev <<std::endl;
00082 
00083    if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent npos2D "<< theEvent->GetNumPos2D() << " nCOPS " << theEvent->GetNumPosCOPS() << std::endl;
00084   
00085   for(int ii=0; ii<theEvent->GetNumPos2D(); ii++) {
00086     AliDaqPosition2D* pos2D = (AliDaqPosition2D*) theEvent->GetArray_Position2D()->At(ii);
00087     if ( ALIUtils::debug >= 4) std::cout<<"2D sensor "<<ii<<" has ID = "<<pos2D->GetID()<< std::endl;
00088       pos2D->DumpIt("2DSENSOR"); 
00089      measList.push_back( GetMeasFromPosition2D( pos2D ) );
00090   }
00091   for(int ii=0; ii<theEvent->GetNumPosCOPS(); ii++) {
00092     AliDaqPositionCOPS* posCOPS = (AliDaqPositionCOPS*) theEvent->GetArray_PositionCOPS()->At(ii);
00093     measList.push_back( GetMeasFromPositionCOPS( posCOPS ) );
00094     if ( ALIUtils::debug >= 4) {
00095       std::cout<<"COPS sensor "<<ii<<" has ID = "<<posCOPS->GetID()<< std::endl;
00096       posCOPS->DumpIt("COPS"); 
00097     }
00098   }
00099   for(int ii=0; ii<theEvent->GetNumTilt(); ii++) {
00100     AliDaqTilt* tilt = (AliDaqTilt*) theEvent->GetArray_Tilt()->At(ii);
00101     measList.push_back( GetMeasFromTilt( tilt ) );
00102      if ( ALIUtils::debug >= 4) {
00103        std::cout<<"TILT sensor "<<ii<<" has ID = "<<tilt->GetID()<< std::endl;
00104        tilt->DumpIt("TILT"); 
00105      }
00106      
00107   }
00108   for(int ii=0; ii<theEvent->GetNumDist(); ii++) {
00109     AliDaqDistance* dist = (AliDaqDistance*) theEvent->GetArray_Dist()->At(ii);
00110     measList.push_back( GetMeasFromDist( dist ) );
00111     if ( ALIUtils::debug >= 4) {
00112       std::cout<<"DIST sensor "<<ii<<" has ID = "<<dist->GetID()<< std::endl;
00113       dist->DumpIt("DIST"); 
00114     }
00115   }
00116   
00117   nextEvent = nev + 1;
00118   
00119   BuildMeasurementsFromOptAlign( measList );
00120   
00121   return 1;
00122   
00123 }
00124 
00125 //----------------------------------------------------------------------
00126 OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromPosition2D( AliDaqPosition2D* pos2D )
00127 {
00128   OpticalAlignMeasurementInfo meas;
00129   
00130   meas.type_ = "SENSOR2D";
00131   meas.name_ = pos2D->GetID();
00132   //-   std::vector<std::string> measObjectNames_;
00133   std::vector<bool> isSimu;
00134   for( size_t jj = 0; jj < 2; jj++ ){
00135     isSimu.push_back(false); 
00136   }
00137   meas.isSimulatedValue_ = isSimu; 
00138   std::vector<OpticalAlignParam> paramList;
00139   OpticalAlignParam oaParam1;
00140   oaParam1.name_ = "H:";
00141   oaParam1.value_ = pos2D->GetX()/100.;
00142   oaParam1.error_ = pos2D->GetXerror()/100.;
00143   paramList.push_back(oaParam1);
00144   
00145   OpticalAlignParam oaParam2;
00146   oaParam2.name_ = "V:";
00147   oaParam2.value_ = pos2D->GetY()/100.;
00148   oaParam2.error_ = pos2D->GetYerror()/100.;
00149   paramList.push_back(oaParam2);
00150   
00151   meas.values_ = paramList;
00152 
00153   return meas;
00154 }
00155 
00156 
00157 //----------------------------------------------------------------------
00158 OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromPositionCOPS( AliDaqPositionCOPS* posCOPS )
00159 {
00160   OpticalAlignMeasurementInfo meas;
00161   
00162   meas.type_ = "COPS";
00163   meas.name_ = posCOPS->GetID();
00164   //-   std::vector<std::string> measObjectNames_;
00165   std::vector<bool> isSimu;
00166   for( size_t jj = 0; jj < 4; jj++ ){
00167     isSimu.push_back(false); 
00168   }
00169   meas.isSimulatedValue_ = isSimu; 
00170 
00171   std::vector<OpticalAlignParam> paramList;
00172   OpticalAlignParam oaParam1;
00173   oaParam1.name_ = "U:";
00174   oaParam1.value_ = posCOPS->GetUp()/100.;
00175   oaParam1.error_ = posCOPS->GetUpError()/100.;
00176   paramList.push_back(oaParam1);
00177 
00178   OpticalAlignParam oaParam2;
00179   oaParam2.name_ = "U:";
00180   oaParam2.value_ = posCOPS->GetDown()/100.;
00181   oaParam2.error_ = posCOPS->GetDownError()/100.;
00182   paramList.push_back(oaParam2);
00183 
00184   OpticalAlignParam oaParam3;
00185   oaParam3.name_ = "U:";
00186   oaParam3.value_ = posCOPS->GetRight()/100.;
00187   oaParam3.error_ = posCOPS->GetRightError()/100.;
00188   paramList.push_back(oaParam3);
00189 
00190   OpticalAlignParam oaParam4;
00191   oaParam4.name_ = "U:";
00192   oaParam4.value_ = posCOPS->GetLeft()/100.;
00193   oaParam4.error_ = posCOPS->GetLeftError()/100.;
00194   paramList.push_back(oaParam4);
00195   
00196   meas.values_ = paramList;
00197 
00198   return meas;
00199 
00200 }
00201 
00202 //----------------------------------------------------------------------
00203 OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromTilt( AliDaqTilt* tilt )
00204 {
00205   OpticalAlignMeasurementInfo meas;
00206   
00207   meas.type_ = "TILTMETER";
00208   meas.name_ = tilt->GetID();
00209   //-   std::vector<std::string> measObjectNames_;
00210   std::vector<bool> isSimu;
00211   for( size_t jj = 0; jj < 2; jj++ ){
00212     isSimu.push_back(false); 
00213   }
00214   meas.isSimulatedValue_ = isSimu; 
00215   std::vector<OpticalAlignParam> paramList;
00216   OpticalAlignParam oaParam;
00217   oaParam.name_ = "T:";
00218   oaParam.value_ = tilt->GetTilt();
00219   oaParam.error_ = tilt->GetTiltError();
00220   paramList.push_back(oaParam);
00221   
00222   meas.values_ = paramList;
00223 
00224   return meas;
00225 
00226 }
00227 
00228 
00229 //----------------------------------------------------------------------
00230 OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromDist( AliDaqDistance* dist )
00231 {
00232   OpticalAlignMeasurementInfo meas;
00233   
00234   meas.type_ = "DISTANCEMETER";
00235   meas.name_ = dist->GetID();
00236   //-   std::vector<std::string> measObjectNames_;
00237   std::vector<bool> isSimu;
00238   for( size_t jj = 0; jj < 2; jj++ ){
00239     isSimu.push_back(false); 
00240   }
00241   meas.isSimulatedValue_ = isSimu; 
00242   std::vector<OpticalAlignParam> paramList;
00243   OpticalAlignParam oaParam;
00244   oaParam.name_ = "D:";
00245   oaParam.value_ = dist->GetDistance()/100.;
00246   oaParam.error_ = dist->GetDistanceError()/100.;
00247   paramList.push_back(oaParam);
00248 
00249   meas.values_ = paramList;
00250 
00251   return meas;
00252 }
00253 
00254 
00255 //----------------------------------------------------------------------
00256 void CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign( std::vector<OpticalAlignMeasurementInfo>& measList )
00257 {
00258   if ( ALIUtils::debug >= 3) std::cout << "@@@ CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign " << std::endl;
00259 
00260  //set date and time of current measurement
00261   //  if( wordlist[0] == "DATE:" ) {
00262   //   Measurement::setCurrentDate( wordlist ); 
00263   // } 
00264 
00265   //---------- loop measurements read from ROOT and check for corresponding measurement in Model
00266   //  ALIint nMeasModel = Model::MeasurementList().size();
00267   ALIint nMeasRoot = measList.size();
00268   if(ALIUtils::debug >= 4) {
00269     std::cout << " Building " << nMeasRoot << " measurements from ROOT file " << std::endl;
00270   }
00271 
00272   //--- Loop to Measurements in Model and check for corresponding measurement in ROOT
00273   std::vector< Measurement* >::const_iterator vmcite;
00274   for( vmcite = Model::MeasurementList().begin();  vmcite != Model::MeasurementList().end(); vmcite++ ) {
00275     ALIint fcolon = (*vmcite)->name().find(':');
00276     ALIstring oname = (*vmcite)->name();
00277     oname = oname.substr(fcolon+1,oname.length());
00278     
00279     //---------- loop measurements read from ROOT 
00280     ALIint ii;
00281     for(ii = 0; ii < nMeasRoot; ii++) {
00282       OpticalAlignMeasurementInfo measInfo = measList[ii];
00283       std::cout << " measurement name ROOT " << measInfo.name_ << " Model= " << (*vmcite)->name() << " short " << oname << std::endl;
00284       
00285       if( oname == measInfo.name_ ) {
00286         //-------- Measurement found, fill data
00287         //---- Check that type is the same
00288         if( (*vmcite)->type() != measInfo.type_ ) {
00289           std::cerr << "!!! Measurement from ROOT file: type in file is " 
00290                     <<measInfo.type_ << " and should be " << (*vmcite)->type() << std::endl;
00291           exit(1);
00292         }
00293         
00294         std::cout << " NOBJECTS IN MEAS " << (*vmcite)->OptOList().size() << " NMEAS " << Model::MeasurementList().size() << std::endl;
00295         
00296         std::vector<OpticalAlignParam> measValues = measInfo.values_;
00297         
00298         for( size_t jj= 0; jj < measValues.size(); jj++ ){
00299           (*vmcite)->fillData( jj, &(measValues[jj]) );
00300         }
00301 
00302         std::cout << " NOBJECTS IN MEAS after " << (*vmcite)->OptOList().size() << " NMEAS " << Model::MeasurementList().size()  << std::endl;
00303 
00304         break;
00305       }
00306     }
00307     if (ii==nMeasRoot) {
00308       std::cerr << "!!! Reading measurement from file: measurement not found! Type in list is "  <<  oname  << std::endl;
00309       exit(1);
00310     }
00311   }
00312 
00313 }
00314