Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301
#include "Alignment/CocoaModel/interface/CocoaDaqReaderRoot.h"
#include "TFile.h"
#include "Alignment/CocoaDaq/interface/CocoaDaqRootEvent.h"
#include "Alignment/CocoaModel/interface/Measurement.h"
#include "Alignment/CocoaModel/interface/Model.h"
#include "Alignment/CocoaUtilities/interface/ALIUtils.h"

#include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurements.h"

#include <iostream>
#include <string>

#include "TClonesArray.h"

//----------------------------------------------------------------------
CocoaDaqReaderRoot::CocoaDaqReaderRoot(const std::string& m_inFileName) {
  if (ALIUtils::debug >= 3)
    std::cout << " CocoaDaqReaderRoot opening file: " << m_inFileName << std::endl;
  // Open root file
  theFile = new TFile(m_inFileName.c_str());
  if (!theTree) {
    std::cerr << " CocoaDaqReaderRoot TTree file not found " << m_inFileName << std::endl;
    throw std::exception();
  }

  // Read TTree named "CocoaDaq" in memory.  !! SHOULD BE CALLED Alignment_Cocoa
  theTree = (TTree*)theFile->Get("CocoaDaq");
  //  theTree = (TTree*)theFile->Get("Alignment_Link_Cocoa");

  if (!theTree) {
    std::cerr << " CocoaDaqReaderRoot TTree in file " << m_inFileName << " should be called 'CocoaDaq' " << std::endl;
    throw std::exception();
  }
  TBranch* branch = theTree->GetBranch("Alignment_Cocoa");

  nev = branch->GetEntries();  // number of entries in Tree
  //if ( ALIUtils::debug >= 2) std::cout << "CocoaDaqReaderRoot::CocoaDaqReaderRoot:  number of entries in Tree " << nev << std::endl;

  nextEvent = 0;

  // Event object must be created before setting the branch address
  theEvent = new CocoaDaqRootEvent();

  // link pointer to Tree branch
  theTree->SetBranchAddress("Alignment_Cocoa", &theEvent);  //  !! SHOULD BE CALLED Alignment_Cocoa
  // theTree->SetBranchAddress("Alignment_Link", &theEvent);  //  !! SHOULD BE CALLED Alignment_Cocoa

  CocoaDaqReader::SetDaqReader(this);
}

//----------------------------------------------------------------------
CocoaDaqReaderRoot::~CocoaDaqReaderRoot() { theFile->Close(); }

//----------------------------------------------------------------------
bool CocoaDaqReaderRoot::ReadNextEvent() { return ReadEvent(nextEvent); }

//----------------------------------------------------------------------
bool CocoaDaqReaderRoot::ReadEvent(int nev) {
  std::vector<OpticalAlignMeasurementInfo> measList;

  int nb = 0;  // dummy, number of bytes
  // Loop over all events
  nb = theTree->GetEntry(nev);  // read in entire event

  if (ALIUtils::debug >= 3)
    std::cout << "CocoaDaqReaderRoot reading event " << nev << " " << nb << std::endl;
  if (nb == 0)
    return false;  //end of file reached??

  // Every n events, dump one to screen
  int n = 1;
  if (nev % n == 0 && ALIUtils::debug >= 3)
    theEvent->DumpIt();

  //if ( ALIUtils::debug >= 3) std::cout<<" CocoaDaqReaderRoot::ReadEvent "<< nev <<std::endl;

  if (ALIUtils::debug >= 3)
    std::cout << " CocoaDaqReaderRoot::ReadEvent npos2D " << theEvent->GetNumPos2D() << " nCOPS "
              << theEvent->GetNumPosCOPS() << std::endl;

  for (int ii = 0; ii < theEvent->GetNumPos2D(); ii++) {
    AliDaqPosition2D* pos2D = (AliDaqPosition2D*)theEvent->GetArray_Position2D()->At(ii);
    if (ALIUtils::debug >= 4)
      std::cout << "2D sensor " << ii << " has ID = " << pos2D->GetID() << std::endl;
    pos2D->DumpIt("2DSENSOR");
    measList.push_back(GetMeasFromPosition2D(pos2D));
  }
  for (int ii = 0; ii < theEvent->GetNumPosCOPS(); ii++) {
    AliDaqPositionCOPS* posCOPS = (AliDaqPositionCOPS*)theEvent->GetArray_PositionCOPS()->At(ii);
    measList.push_back(GetMeasFromPositionCOPS(posCOPS));
    if (ALIUtils::debug >= 4) {
      std::cout << "COPS sensor " << ii << " has ID = " << posCOPS->GetID() << std::endl;
      posCOPS->DumpIt("COPS");
    }
  }
  for (int ii = 0; ii < theEvent->GetNumTilt(); ii++) {
    AliDaqTilt* tilt = (AliDaqTilt*)theEvent->GetArray_Tilt()->At(ii);
    measList.push_back(GetMeasFromTilt(tilt));
    if (ALIUtils::debug >= 4) {
      std::cout << "TILT sensor " << ii << " has ID = " << tilt->GetID() << std::endl;
      tilt->DumpIt("TILT");
    }
  }
  for (int ii = 0; ii < theEvent->GetNumDist(); ii++) {
    AliDaqDistance* dist = (AliDaqDistance*)theEvent->GetArray_Dist()->At(ii);
    measList.push_back(GetMeasFromDist(dist));
    if (ALIUtils::debug >= 4) {
      std::cout << "DIST sensor " << ii << " has ID = " << dist->GetID() << std::endl;
      dist->DumpIt("DIST");
    }
  }

  nextEvent = nev + 1;

  BuildMeasurementsFromOptAlign(measList);

  return true;
}

//----------------------------------------------------------------------
OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromPosition2D(AliDaqPosition2D* pos2D) {
  OpticalAlignMeasurementInfo meas;

  meas.type_ = "SENSOR2D";
  meas.name_ = std::string(pos2D->GetID().Data());
  //-   std::vector<std::string> measObjectNames_;
  std::vector<bool> isSimu;
  for (size_t jj = 0; jj < 2; jj++) {
    isSimu.push_back(false);
  }
  meas.isSimulatedValue_ = isSimu;
  std::vector<OpticalAlignParam> paramList;
  OpticalAlignParam oaParam1;
  oaParam1.name_ = "H:";
  oaParam1.value_ = pos2D->GetX() / 100.;
  oaParam1.error_ = pos2D->GetXerror() / 100.;
  paramList.push_back(oaParam1);

  OpticalAlignParam oaParam2;
  oaParam2.name_ = "V:";
  oaParam2.value_ = pos2D->GetY() / 100.;
  oaParam2.error_ = pos2D->GetYerror() / 100.;
  paramList.push_back(oaParam2);

  meas.values_ = paramList;

  return meas;
}

//----------------------------------------------------------------------
OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromPositionCOPS(AliDaqPositionCOPS* posCOPS) {
  OpticalAlignMeasurementInfo meas;

  meas.type_ = "COPS";
  meas.name_ = std::string(posCOPS->GetID().Data());
  //-   std::vector<std::string> measObjectNames_;
  std::vector<bool> isSimu;
  for (size_t jj = 0; jj < 4; jj++) {
    isSimu.push_back(false);
  }
  meas.isSimulatedValue_ = isSimu;

  std::vector<OpticalAlignParam> paramList;
  OpticalAlignParam oaParam1;
  oaParam1.name_ = "U:";
  oaParam1.value_ = posCOPS->GetUp() / 100.;
  oaParam1.error_ = posCOPS->GetUpError() / 100.;
  paramList.push_back(oaParam1);

  OpticalAlignParam oaParam2;
  oaParam2.name_ = "U:";
  oaParam2.value_ = posCOPS->GetDown() / 100.;
  oaParam2.error_ = posCOPS->GetDownError() / 100.;
  paramList.push_back(oaParam2);

  OpticalAlignParam oaParam3;
  oaParam3.name_ = "U:";
  oaParam3.value_ = posCOPS->GetRight() / 100.;
  oaParam3.error_ = posCOPS->GetRightError() / 100.;
  paramList.push_back(oaParam3);

  OpticalAlignParam oaParam4;
  oaParam4.name_ = "U:";
  oaParam4.value_ = posCOPS->GetLeft() / 100.;
  oaParam4.error_ = posCOPS->GetLeftError() / 100.;
  paramList.push_back(oaParam4);

  meas.values_ = paramList;

  return meas;
}

//----------------------------------------------------------------------
OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromTilt(AliDaqTilt* tilt) {
  OpticalAlignMeasurementInfo meas;

  meas.type_ = "TILTMETER";
  meas.name_ = std::string(tilt->GetID().Data());
  //-   std::vector<std::string> measObjectNames_;
  std::vector<bool> isSimu;
  for (size_t jj = 0; jj < 2; jj++) {
    isSimu.push_back(false);
  }
  meas.isSimulatedValue_ = isSimu;
  std::vector<OpticalAlignParam> paramList;
  OpticalAlignParam oaParam;
  oaParam.name_ = "T:";
  oaParam.value_ = tilt->GetTilt();
  oaParam.error_ = tilt->GetTiltError();
  paramList.push_back(oaParam);

  meas.values_ = paramList;

  return meas;
}

//----------------------------------------------------------------------
OpticalAlignMeasurementInfo CocoaDaqReaderRoot::GetMeasFromDist(AliDaqDistance* dist) {
  OpticalAlignMeasurementInfo meas;

  meas.type_ = "DISTANCEMETER";
  meas.name_ = std::string(dist->GetID().Data());
  //-   std::vector<std::string> measObjectNames_;
  std::vector<bool> isSimu;
  for (size_t jj = 0; jj < 2; jj++) {
    isSimu.push_back(false);
  }
  meas.isSimulatedValue_ = isSimu;
  std::vector<OpticalAlignParam> paramList;
  OpticalAlignParam oaParam;
  oaParam.name_ = "D:";
  oaParam.value_ = dist->GetDistance() / 100.;
  oaParam.error_ = dist->GetDistanceError() / 100.;
  paramList.push_back(oaParam);

  meas.values_ = paramList;

  return meas;
}

//----------------------------------------------------------------------
void CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign(std::vector<OpticalAlignMeasurementInfo>& measList) {
  if (ALIUtils::debug >= 3)
    std::cout << "@@@ CocoaDaqReaderRoot::BuildMeasurementsFromOptAlign " << std::endl;

  //set date and time of current measurement
  //  if( wordlist[0] == "DATE:" ) {
  //   Measurement::setCurrentDate( wordlist );
  // }

  //---------- loop measurements read from ROOT and check for corresponding measurement in Model
  //  ALIint nMeasModel = Model::MeasurementList().size();
  ALIint nMeasRoot = measList.size();
  if (ALIUtils::debug >= 4) {
    std::cout << " Building " << nMeasRoot << " measurements from ROOT file " << std::endl;
  }

  //--- Loop to Measurements in Model and check for corresponding measurement in ROOT
  std::vector<Measurement*>::const_iterator vmcite;
  for (vmcite = Model::MeasurementList().begin(); vmcite != Model::MeasurementList().end(); ++vmcite) {
    ALIint fcolon = (*vmcite)->name().find(':');
    ALIstring oname = (*vmcite)->name();
    oname = oname.substr(fcolon + 1, oname.length());

    //---------- loop measurements read from ROOT
    ALIint ii;
    for (ii = 0; ii < nMeasRoot; ii++) {
      OpticalAlignMeasurementInfo measInfo = measList[ii];
      std::cout << " measurement name ROOT " << measInfo.name_ << " Model= " << (*vmcite)->name() << " short " << oname
                << std::endl;

      if (oname == measInfo.name_) {
        //-------- Measurement found, fill data
        //---- Check that type is the same
        if ((*vmcite)->type() != measInfo.type_) {
          std::cerr << "!!! Measurement from ROOT file: type in file is " << measInfo.type_ << " and should be "
                    << (*vmcite)->type() << std::endl;
          exit(1);
        }

        std::cout << " NOBJECTS IN MEAS " << (*vmcite)->OptOList().size() << " NMEAS "
                  << Model::MeasurementList().size() << std::endl;

        std::vector<OpticalAlignParam> measValues = measInfo.values_;

        for (size_t jj = 0; jj < measValues.size(); jj++) {
          (*vmcite)->fillData(jj, &(measValues[jj]));
        }

        std::cout << " NOBJECTS IN MEAS after " << (*vmcite)->OptOList().size() << " NMEAS "
                  << Model::MeasurementList().size() << std::endl;

        break;
      }
    }
    if (ii == nMeasRoot) {
      std::cerr << "!!! Reading measurement from file: measurement not found! Type in list is " << oname << std::endl;
      exit(1);
    }
  }
}