FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndLmdAlignQA Class Reference

#include <PndLmdAlignQA.h>

Public Member Functions

 PndLmdAlignQA ()
 
virtual ~PndLmdAlignQA ()
 
void init ()
 
void checkCombined ()
 
void checkCyclicMatrices (bool inCentimeters=true)
 
void checkCombinedMatrices (bool inCentimeters=true)
 
std::map< std::string,
TGeoHMatrix > * 
readRootMatrices (TString &filename)
 
TGeoHMatrix baseTransformation (TGeoHMatrix &input, TGeoHMatrix &toBaseMatrix)
 
TGeoHMatrix getMatrixSensorToSensor (int sensorOne, int sensorTwo)
 
TGeoHMatrix getOverlapMatrixLikeICP (PndLmdOverlapInfo &info)
 
TGeoHMatrix getMisalignedOverlapFromGeoManager (PndLmdOverlapInfo &info)
 
TGeoHMatrix getMisalignedOverlapFromICP (PndLmdOverlapInfo &info, std::string ICPmatrix)
 
std::vector< double > getMatrixDiffCM (PndLmdOverlapInfo &info, std::string &icpFile)
 
PndLmdOverlapInfogetSmallOverlapInfo (std::vector< PndLmdOverlapInfo > &infos, int smallOverlap)
 
void readMatrixInfo ()
 
void checkIOpaths ()
 
bool checkForMatrixFiles ()
 
void calculateOverlapingAreas ()
 
void setInCentimeters (bool inCentimeters)
 
void setLmdMatPath (const std::string &path)
 
void setBinaryMatPath (const std::string &path)
 
void setPdfOutPath (const std::string &path)
 
void setAlignedGeometry (bool aligned)
 
void setPairsRequired (int number)
 
void setOutputPath (std::string path)
 
void setInfoAbsolute (bool info)
 
void setInfoMomentum (double info)
 
void setInfoRelative (bool info)
 

Private Member Functions

void createHist (std::vector< std::vector< double > > &vec, histParams &parameters)
 
double calculateOverlappingArea (int id1, int id2)
 
int noOfPairs (int overlapID)
 

Private Attributes

std::vector< std::string > _inputFiles
 
std::string outputPath
 
std::string pdfOutPath
 
std::string binaryMatPath
 
std::string LMDMatPath
 
double infoMomentum
 
bool infoAbsolute
 
bool infoRelative
 
bool byPlane
 
bool _inCentimeters
 
bool alignOptionBool
 
int pairsRequired
 
PndLmdAlignManager manager
 
PndLmdGeometryHelperhelper
 
std::map< int, int > matrixInfo
 
std::map< std::string,
TGeoHMatrix > * 
matricesMisaligned
 

Detailed Description

Definition at line 34 of file PndLmdAlignQA.h.

Constructor & Destructor Documentation

PndLmdAlignQA::PndLmdAlignQA ( )

Definition at line 36 of file PndLmdAlignQA.cxx.

References init().

36  {
37  init();
38 }
PndLmdAlignQA::~PndLmdAlignQA ( )
virtual

Definition at line 40 of file PndLmdAlignQA.cxx.

40  {
41 
42 }

Member Function Documentation

TGeoHMatrix PndLmdAlignQA::baseTransformation ( TGeoHMatrix &  input,
TGeoHMatrix &  toBaseMatrix 
)

Definition at line 154 of file PndLmdAlignQA.cxx.

Referenced by getOverlapMatrixLikeICP().

154  {
155  return TGeoHMatrix(toBaseMatrix * input * toBaseMatrix.Inverse());
156 }
void PndLmdAlignQA::calculateOverlapingAreas ( )

Definition at line 245 of file PndLmdAlignQA.cxx.

References calculateOverlappingArea(), exit(), PndLmdGeometryHelper::getAvailableOverlapIDs(), PndLmdGeometryHelper::getSensorOneFromOverlapID(), PndLmdGeometryHelper::getSensorTwoFromOverlapID(), helper, and i.

245  {
246 
247  //get all available overlapping areas
248 
249  cout << "------------ TESTING --------------\n";
250  cout << "\\AtoB{4}{9} & " << calculateOverlappingArea(4, 9) << "\n";
251  cout << "\\AtoB{0}{1} & " << calculateOverlappingArea(0, 1) << "\n";
252  cout << "\\AtoB{0}{2} & " << calculateOverlappingArea(0, 2) << "\n";
253  cout << "\\AtoB{0}{9} & " << calculateOverlappingArea(0, 9) << "\n";
254  cout << "------------ DONE --------------\n";
255 
256  exit(0);
257 
258  vector<int> overlapIDs = helper->getAvailableOverlapIDs();
259  int id1, id2;
260  double areaPercent = 0;
261 
262  //get id1 and id2 from them and calc
263  for (unsigned int i = 0; i < overlapIDs.size(); i++) {
264  id1 = helper->getSensorOneFromOverlapID(overlapIDs[i]);
265  id2 = helper->getSensorTwoFromOverlapID(overlapIDs[i]);
266 
267  areaPercent = calculateOverlappingArea(id1, id2);
268  cout << "\\AtoB{" << id1 << "}{" << id2 << "} & " << areaPercent << "\n";
269  }
270  exit(0);
271 
272 }
std::vector< int > getAvailableOverlapIDs()
int getSensorTwoFromOverlapID(int overlapID)
Int_t i
Definition: run_full.C:25
exit(0)
double calculateOverlappingArea(int id1, int id2)
int getSensorOneFromOverlapID(int overlapID)
PndLmdGeometryHelper * helper
Definition: PndLmdAlignQA.h:45
double PndLmdAlignQA::calculateOverlappingArea ( int  id1,
int  id2 
)
private

Definition at line 204 of file PndLmdAlignQA.cxx.

References getMatrixSensorToSensor().

Referenced by calculateOverlapingAreas().

204  {
205 
206  TGeoHMatrix sen1to2 = getMatrixSensorToSensor(sensor1, sensor2);
207 
208  //histogram those distances
209  double colTest2 = 0, rowTest2 = 0;
210  int valid = 0;
211 
212  TVector3 testVector;
213  TVector3 testVectorTransformed;
214 
215  //for every pixel (250*250) find all 4 overlap pixels (cutoff at 100 microns, pixels can't be that far away on perfect geometry)
216  for (int colTest1 = 0; colTest1 < 247; colTest1++) {
217  for (int rowTest1 = 0; rowTest1 < 242; rowTest1++) {
218 
219  testVector.SetXYZ(colTest1, rowTest1, 0.0);
220 
221  //FIXME: this line is not correct
222  //sen1to2.MasterToLocal(testVector., testVectorTransformed);
223 
224  colTest2 = testVectorTransformed.x();
225  rowTest2 = testVectorTransformed.y();
226 
227  //also, check if overlap pixel even exists. must be row elem [0,250], col elem [0,250]
228  //are we still overlapping area?
229  if (colTest2 < 0 || colTest2 > 247) { //to account for inactive area
230  continue;
231  }
232  if (rowTest2 < 0 || rowTest2 > 242) { //to account for inactive area
233  continue;
234  }
235 
236  //we are still on overlapping area
237  valid++;
238  }
239  }
240  double coverage = (double) valid / (250.0 * 250.0) * 100;
241 
242  return coverage;
243 }
TGeoHMatrix getMatrixSensorToSensor(int sensorOne, int sensorTwo)
TGeoVolume * sensor2
void PndLmdAlignQA::checkCombined ( )

Definition at line 370 of file PndLmdAlignQA.cxx.

References histParams::bins, createHist(), ext, histParams::fileName, getMatrixDiff(), getMatrixSensorToSensor(), getMisalignedOverlapFromGeoManager(), getMisalignedOverlapFromICP(), PndLmdGeometryHelper::getOverlapInfos(), getSmallOverlapInfo(), helper, i, PndLmdOverlapInfo::id1, PndLmdOverlapInfo::id2, manager, matricesMisaligned, histParams::path, pdfOutPath, readRootMatrices(), histParams::scaleFactor, sensor2, PndLmdAlignManager::setMatrixOutDir(), histParams::title, TString, histParams::vectorIndex, histParams::xtitle, and histParams::ytitle.

370  {
371 
372  //TODO: make more general, set paths during setup
373  TString misalignedMatrices = "misalignMatrices-SensorsOnly-100u.root";
374  TString idealMatrices = "idealMatrices.root";
375  std::string path =
376  "/home/arbeit/RedPro3TB/simulationData/2018-05-07-misalign-100u/LMDmatrices-inSensorOne";
377 
378  matricesMisaligned = readRootMatrices(misalignedMatrices);
379 
380  std::map<int, TGeoHMatrix> matricesByHand;
381  std::map<int, TGeoHMatrix> matricesByGeoManager;
382 
383  bool compareWithIdeal = true;
384  string ext = "cm.mat";
385 
386  cout << "getting matrices by hand/ICP-like\n";
387  for (int iHalf = 0; iHalf < 2; iHalf++) {
388  for (int iPlane = 0; iPlane < 4; iPlane++) {
389  for (int iModule = 0; iModule < 5; iModule++) {
390 
391  // the order in the vector is not defined, use getSmallOverlapInfo() for that
392  std::vector<PndLmdOverlapInfo> overlaps = helper->getOverlapInfos(iHalf, iPlane, iModule);
393 
394  // load all 9 overlap matrices
395  TGeoHMatrix mat0;
396  TGeoHMatrix mat1;
397  TGeoHMatrix mat2;
398  TGeoHMatrix mat3;
399  TGeoHMatrix mat4;
400  TGeoHMatrix mat5;
401  TGeoHMatrix mat6;
402  TGeoHMatrix mat7;
403  TGeoHMatrix mat8;
404 
405  if (compareWithIdeal) {
415  }
416  else {
417  string icp0 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 0).overlapID) + ext;
418  mat0 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 0), icp0);
419  string icp1 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 1).overlapID) + ext;
420  mat1 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 1), icp1);
421  string icp2 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 2).overlapID) + ext;
422  mat2 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 2), icp2);
423  string icp3 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 3).overlapID) + ext;
424  mat3 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 3), icp3);
425  string icp4 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 4).overlapID) + ext;
426  mat4 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 4), icp4);
427  string icp5 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 5).overlapID) + ext;
428  mat5 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 5), icp5);
429  string icp6 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 6).overlapID) + ext;
430  mat6 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 6), icp6);
431  string icp7 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 7).overlapID) + ext;
432  mat7 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 7), icp7);
433  string icp8 = path + "/m" + std::to_string(getSmallOverlapInfo(overlaps, 8).overlapID) + ext;
434  mat8 = getMisalignedOverlapFromICP(getSmallOverlapInfo(overlaps, 8), icp8);
435  }
436 
437  // get 0->1 and 5->6
438  TGeoHMatrix sen0to1 = getMatrixSensorToSensor(0, 1);
439  TGeoHMatrix sen5to6 = getMatrixSensorToSensor(5, 6);
440 
441  TGeoHMatrix sen1to2 = mat4 * mat5.Inverse();
442  TGeoHMatrix sen1to3 = mat4 * mat1.Inverse();
443  TGeoHMatrix sen1to4a = mat4 * mat5.Inverse() * mat6 * mat2.Inverse();
444  TGeoHMatrix sen1to4b = mat4 * mat1.Inverse() * mat7 * mat8.Inverse();
445 
446  TGeoHMatrix sen1to5 = mat0; // well, kind of...
447  TGeoHMatrix sen1to6 = mat4 * mat1.Inverse() * mat3;
448  TGeoHMatrix sen1to7 = mat4 * mat1.Inverse() * mat7;
449  TGeoHMatrix sen1to8 = mat4;
450  TGeoHMatrix sen1to9a = mat4 * mat5.Inverse() * mat6;
451  TGeoHMatrix sen1to9b = mat4 * mat1.Inverse() * mat7 * mat8.Inverse() * mat2;
452 
453  int matrixId = 1000 * iHalf + 100 * iPlane + 10 * iModule;
454 
455  matricesByHand[matrixId + 2] = sen1to2;
456  matricesByHand[matrixId + 3] = sen1to3;
457  matricesByHand[matrixId + 4] = sen1to4a;
458  matricesByHand[matrixId + 5] = sen1to5; // this is not yet correct
459  matricesByHand[matrixId + 6] = sen1to6;
460  matricesByHand[matrixId + 7] = sen1to7;
461  matricesByHand[matrixId + 8] = sen1to8;
462  matricesByHand[matrixId + 9] = sen1to9a;
463  }
464  }
465  }
466 
467  // TODO: this must be implemented
468  //applyMisalignmentToGGeoManager();
469  cout << "getting matrices by geoManager after misalignment\n";
470  for (int iHalf = 0; iHalf < 2; iHalf++) {
471  for (int iPlane = 0; iPlane < 4; iPlane++) {
472  for (int iModule = 0; iModule < 5; iModule++) {
473 
474  // the order in the vector is not defined, use getSmallOverlapInfo() for that
475  std::vector<PndLmdOverlapInfo> overlaps = helper->getOverlapInfos(iHalf, iPlane, iModule);
476  // get everything from gGeoManager
477  int matrixId = 1000 * iHalf + 100 * iPlane + 10 * iModule;
478 
479  PndLmdOverlapInfo info0 = getSmallOverlapInfo(overlaps, 0);
480  PndLmdOverlapInfo info1 = getSmallOverlapInfo(overlaps, 1);
481  PndLmdOverlapInfo info2 = getSmallOverlapInfo(overlaps, 2);
482  PndLmdOverlapInfo info3 = getSmallOverlapInfo(overlaps, 3);
483  PndLmdOverlapInfo info4 = getSmallOverlapInfo(overlaps, 4);
484  PndLmdOverlapInfo info5 = getSmallOverlapInfo(overlaps, 5);
485  PndLmdOverlapInfo info6 = getSmallOverlapInfo(overlaps, 6);
486  PndLmdOverlapInfo info7 = getSmallOverlapInfo(overlaps, 7);
487  PndLmdOverlapInfo info8 = getSmallOverlapInfo(overlaps, 8);
488 
489  // get sensors 1 and 2 from the overlapInfos
490  int sensor0 = info0.id1;
491  int sensor1 = info4.id1;
492  int sensor2 = info5.id1;
493  int sensor3 = info1.id1;
494  int sensor4 = info2.id1;
495  int sensor5 = info0.id2;
496  int sensor6 = info3.id2;
497  int sensor7 = info7.id2;
498  int sensor8 = info4.id2;
499  int sensor9 = info6.id2;
500 
501  TGeoHMatrix sen1to2 = getMatrixSensorToSensor(sensor1, sensor2);
502  TGeoHMatrix sen1to3 = getMatrixSensorToSensor(sensor1, sensor3);
503  TGeoHMatrix sen1to4a = getMatrixSensorToSensor(sensor1, sensor4);
504  TGeoHMatrix sen1to5 = getMatrixSensorToSensor(sensor0, sensor5);
505  TGeoHMatrix sen1to6 = getMatrixSensorToSensor(sensor1, sensor6);
506  TGeoHMatrix sen1to7 = getMatrixSensorToSensor(sensor1, sensor7);
507  TGeoHMatrix sen1to8 = getMatrixSensorToSensor(sensor1, sensor8);
508  TGeoHMatrix sen1to9a = getMatrixSensorToSensor(sensor1, sensor9);
509 
510  matricesByGeoManager[matrixId + 2] = sen1to2;
511  matricesByGeoManager[matrixId + 3] = sen1to3;
512  matricesByGeoManager[matrixId + 4] = sen1to4a;
513  matricesByGeoManager[matrixId + 5] = sen1to5; // this is not yet correct
514  matricesByGeoManager[matrixId + 6] = sen1to6;
515  matricesByGeoManager[matrixId + 7] = sen1to7;
516  matricesByGeoManager[matrixId + 8] = sen1to8;
517  matricesByGeoManager[matrixId + 9] = sen1to9a;
518  }
519  }
520  }
521 
522  cout << "making histograms\n";
523 
524  //data holds sets of entries. order is (matrixID, alpha, dx, dy
525  std::vector<std::vector<double> > data;
526 
527  for (auto &i : matricesByHand) {
528  int key = i.first;
529  auto thisDiff = getMatrixDiff(matricesByHand[key], matricesByGeoManager[key]);
530 
531  std::vector<double> result;
532  result.push_back(key);
533  result.push_back(thisDiff[2]); // sin(alpha)
534  result.push_back(thisDiff[0]); // tx
535  result.push_back(thisDiff[1]); // ty
536  data.push_back(result);
537  }
538 
539  //prepare
540  string pdfdir = pdfOutPath;
541  manager.setMatrixOutDir(path);
542 
543  histParams parameters;
544  parameters.bins = 20;
545  parameters.path = pdfdir + "/residualsCombined/";
546  parameters.ytitle = "entries";
547 
548  //for DX
549  parameters.title = "matrix combined 0-9 CM - Target, #DeltaX";
550  parameters.xtitle = "dX [#mum]";
551  parameters.scaleFactor = 1e4;
552  parameters.fileName = "dx.pdf";
553  parameters.vectorIndex = 2;
554  createHist(data, parameters);
555 
556  //for DY
557  parameters.title = "matrix combined 0-9 CM - Target, #DeltaY";
558  parameters.xtitle = "dY [#mum]";
559  parameters.scaleFactor = 1e4;
560  parameters.fileName = "dy.pdf";
561  parameters.vectorIndex = 3;
562  createHist(data, parameters);
563 
564  //for DAlpha
565  parameters.title = "matrix combined 0-9 CM - Target, #Delta#alpha";
566  parameters.xtitle = "d#alpha [#murad]";
567  parameters.scaleFactor = 1e6;
568  parameters.fileName = "dalpha.pdf";
569  parameters.vectorIndex = 1;
570  createHist(data, parameters);
571 
572 }
void setMatrixOutDir(std::string directory)
std::map< std::string, TGeoHMatrix > * matricesMisaligned
Definition: PndLmdAlignQA.h:51
PndLmdOverlapInfo & getSmallOverlapInfo(std::vector< PndLmdOverlapInfo > &infos, int smallOverlap)
Int_t i
Definition: run_full.C:25
TGeoHMatrix getMisalignedOverlapFromICP(PndLmdOverlapInfo &info, std::string ICPmatrix)
std::string title
Definition: PndLmdAlignQA.h:27
std::vector< double > getMatrixDiff(TGeoHMatrix &mat1, TGeoHMatrix &mat2)
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
double scaleFactor
Definition: PndLmdAlignQA.h:26
std::vector< PndLmdOverlapInfo > getOverlapInfos(int iHalf=-1, int iPlane=-1, int iModule=-1)
std::string pdfOutPath
Definition: PndLmdAlignQA.h:39
std::string ytitle
Definition: PndLmdAlignQA.h:27
TGeoHMatrix getMatrixSensorToSensor(int sensorOne, int sensorTwo)
PndLmdGeometryHelper * helper
Definition: PndLmdAlignQA.h:45
TGeoVolume * sensor2
std::string path
Definition: PndLmdAlignQA.h:27
TNtuple * ext
Definition: reco_muo.C:24
std::string xtitle
Definition: PndLmdAlignQA.h:27
void createHist(std::vector< std::vector< double > > &vec, histParams &parameters)
TGeoHMatrix getMisalignedOverlapFromGeoManager(PndLmdOverlapInfo &info)
std::map< std::string, TGeoHMatrix > * readRootMatrices(TString &filename)
std::string fileName
Definition: PndLmdAlignQA.h:27
void PndLmdAlignQA::checkCombinedMatrices ( bool  inCentimeters = true)
void PndLmdAlignQA::checkCyclicMatrices ( bool  inCentimeters = true)

Definition at line 576 of file PndLmdAlignQA.cxx.

References _inCentimeters, histParams::bins, createHist(), histParams::fileName, i, LMDMatPath, manager, histParams::path, pdfOutPath, histParams::printCMPXinPathName, readMatrixInfo(), histParams::scaleFactor, PndLmdAlignManager::setInCentimeters(), PndLmdAlignManager::setMatrixOutDir(), histParams::title, Matrix::val, histParams::vectorIndex, histParams::xMax, histParams::xMin, histParams::xtitle, and histParams::ytitle.

576  {
577 
578  if (inCentimeters) {
579  _inCentimeters = true;
581  }
582  else {
583  _inCentimeters = false;
585  }
586 
587  //data holds sets of entries. order is (id1, id2, alpha, dx, dy
588  std::vector<std::vector<double> > data;
589  readMatrixInfo();
590 
591  //prepare
592  string path = LMDMatPath;
593  string pdfdir = pdfOutPath;
594  manager.setMatrixOutDir(path);
595 
596  Matrix cycle;
597  for (int i = 0; i < 400; i++) {
598 
599  // does not matter if in LMC local or sensor local, should always be identity matrix!
600  //cycle = manager.combineCyclicMatrix(i);
601 
602  //store this residual tuple to data
603  std::vector<double> result;
604  result.push_back(0);
605  result.push_back(0);
606  result.push_back(cycle.val[0][1]); // sin(alpha)
607  result.push_back(cycle.val[0][3]); // tx
608  result.push_back(cycle.val[1][3]); // ty
609  data.push_back(result);
610  }
611 
612  histParams parameters;
613 
614  parameters.bins = 20;
615 
616  //for DX
617  parameters.path = pdfdir + "/cyclicChecks/";
618  inCentimeters ? parameters.title = "matrixCM cycle check, #DeltaX" : parameters.title =
619  "matrixPX cycle check, #DeltaX";
620  parameters.xtitle = "dX [#mum]";
621  parameters.ytitle = "entries";
622  parameters.scaleFactor = 1e4;
623  parameters.fileName = "dx.pdf";
624  parameters.vectorIndex = 3;
625  parameters.xMin = -1;
626  parameters.xMax = -1;
627  parameters.printCMPXinPathName = true;
628  createHist(data, parameters);
629 
630  //for DY
631  //parameters.path = pdfdir;
632  inCentimeters ? parameters.title = "matrixCM cycle check, #DeltaY" : parameters.title =
633  "matrixPX cycle check, #DeltaY";
634  parameters.xtitle = "dY [#mum]";
635  parameters.ytitle = "entries";
636  parameters.scaleFactor = 1e4;
637  parameters.fileName = "dy.pdf";
638  parameters.vectorIndex = 4;
639  parameters.xMin = -1;
640  parameters.xMax = -1;
641  parameters.printCMPXinPathName = true;
642  createHist(data, parameters);
643 
644  //for DAlpha
645  //parameters.path = pdfdir;
646  inCentimeters ? parameters.title = "matrixCM cycle check, #Delta#alpha" : parameters.title =
647  "matrixPX cycle check, #Delta#alpha";
648  parameters.xtitle = "d#alpha [#murad]";
649  parameters.ytitle = "entries";
650  parameters.scaleFactor = 1e6;
651  parameters.fileName = "dalpha.pdf";
652  parameters.vectorIndex = 2;
653  parameters.xMin = -1;
654  parameters.xMax = -1;
655  parameters.printCMPXinPathName = true;
656  createHist(data, parameters);
657 }
void setMatrixOutDir(std::string directory)
void setInCentimeters(bool inCentimeters)
Int_t i
Definition: run_full.C:25
std::string title
Definition: PndLmdAlignQA.h:27
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
std::string LMDMatPath
Definition: PndLmdAlignQA.h:39
double scaleFactor
Definition: PndLmdAlignQA.h:26
double xMin
Definition: PndLmdAlignQA.h:29
std::string pdfOutPath
Definition: PndLmdAlignQA.h:39
double xMax
Definition: PndLmdAlignQA.h:30
Definition: matrix.h:50
std::string ytitle
Definition: PndLmdAlignQA.h:27
void readMatrixInfo()
std::string path
Definition: PndLmdAlignQA.h:27
std::string xtitle
Definition: PndLmdAlignQA.h:27
void createHist(std::vector< std::vector< double > > &vec, histParams &parameters)
bool printCMPXinPathName
Definition: PndLmdAlignQA.h:25
std::string fileName
Definition: PndLmdAlignQA.h:27
FLOAT ** val
Definition: matrix.h:138
bool PndLmdAlignQA::checkForMatrixFiles ( )

Definition at line 278 of file PndLmdAlignQA.cxx.

References _inCentimeters, files, PndLmdGeometryHelper::getAvailableOverlapIDs(), helper, i, LMDMatPath, PndLmdAlignManager::makeMatrixFileName(), manager, and PndLmdAlignManager::searchFiles().

Referenced by runLumiPixel2gAlignQA().

278  {
279 
280  //list all IDs that SHOULD be there
281  vector<int> availableIds = helper->getAvailableOverlapIDs();
282 
283  vector<string> files;
284  manager.searchFiles(LMDMatPath, files, "mat", false);
285  int foundFiles = 0;
286 
287  //no matrix files at all!
288  if (files.size() == 0) {
289  return false;
290  }
291 
292  string matrixName;
293  bool tempfilefound = false;
294 
295  //check for every ID that should be there if there is a corresponding file
296  for (size_t i = 0; i < availableIds.size(); i++) {
297 
298  //reset counter
299  tempfilefound = false;
300  matrixName = manager.makeMatrixFileName(availableIds[i], _inCentimeters);
301 
302  for (size_t j = 0; j < files.size(); j++) {
303  if (files[j].find(matrixName) != string::npos) {
304  tempfilefound = true;
305  foundFiles++;
306  }
307  }
308  //file not found? at least one is missing, return false
309  if (!tempfilefound) {
310  return tempfilefound;
311  }
312  }
313  // if no file could not be found, everything is okay
314  return true;
315 }
std::vector< int > getAvailableOverlapIDs()
Int_t i
Definition: run_full.C:25
static int searchFiles(std::string curr_directory, std::vector< std::string > &list, std::string extension="", bool includeSubDirs=true)
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
std::map< int, TString > files
Definition: simubg.C:28
std::string LMDMatPath
Definition: PndLmdAlignQA.h:39
PndLmdGeometryHelper * helper
Definition: PndLmdAlignQA.h:45
static std::string makeMatrixFileName(int overlapId=0, bool incentimeters=true)
void PndLmdAlignQA::checkIOpaths ( )
void PndLmdAlignQA::createHist ( std::vector< std::vector< double > > &  vec,
histParams parameters 
)
private

Definition at line 57 of file PndLmdAlignQA.cxx.

References _inCentimeters, histParams::bins, exit(), histParams::fileName, PndLmdAlignManager::mkdir(), histParams::path, histParams::printCMPXinPathName, histParams::scaleFactor, histParams::title, vec, histParams::vectorIndex, histParams::xMax, histParams::xMin, histParams::xtitle, and histParams::ytitle.

Referenced by checkCombined(), and checkCyclicMatrices().

57  {
58 
59  TH1D histogram(parameters.title.c_str(), parameters.title.c_str(), parameters.bins, parameters.xMin,
60  parameters.xMax);
61 
62  histogram.GetXaxis()->SetTitle(parameters.xtitle.c_str());
63  histogram.GetYaxis()->SetTitle(parameters.ytitle.c_str());
64  if (vec.size() == 0) {
65  cout << "Error: nothing read, data empty! (maybe not enough pairs?) \n";
66  exit(1);
67  }
68 
69  double dataPoint;
70  for (size_t iArea = 0; iArea < vec.size(); iArea++) {
71  dataPoint = vec[iArea][parameters.vectorIndex] * parameters.scaleFactor;
72  histogram.Fill(dataPoint);
73  }
74  std::stringstream pathname;
75  pathname << parameters.path;
76 
77  if (parameters.printCMPXinPathName) {
78  _inCentimeters ? pathname << "/inCm/" : pathname << "/inPx/";
79  }
80  PndLmdAlignManager::mkdir(pathname.str());
81  TCanvas canvas("canvas", "canvas", 800, 600);
82  canvas.cd();
83  histogram.Draw();
84  histogram.Draw("HIST TEXT0 SAME");
85  canvas.Print((pathname.str() + parameters.fileName).c_str());
86 }
exit(0)
std::string title
Definition: PndLmdAlignQA.h:27
double scaleFactor
Definition: PndLmdAlignQA.h:26
static bool mkdir(std::string path)
double xMin
Definition: PndLmdAlignQA.h:29
double xMax
Definition: PndLmdAlignQA.h:30
std::string ytitle
Definition: PndLmdAlignQA.h:27
std::string path
Definition: PndLmdAlignQA.h:27
std::string xtitle
Definition: PndLmdAlignQA.h:27
bool printCMPXinPathName
Definition: PndLmdAlignQA.h:25
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
std::string fileName
Definition: PndLmdAlignQA.h:27
std::vector< double > PndLmdAlignQA::getMatrixDiffCM ( PndLmdOverlapInfo info,
std::string &  icpFile 
)

Definition at line 317 of file PndLmdAlignQA.cxx.

References dx, dy, getOverlapMatrixLikeICP(), manager, and PndLmdAlignManager::readTGeoHMatrix().

317  {
318  std::vector<double> result;
319 
320  TGeoHMatrix ICP = manager.readTGeoHMatrix(icpFile);
321  TGeoHMatrix Mtarget = getOverlapMatrixLikeICP(info);
322 
323  const double *translationDesign = Mtarget.GetTranslation();
324  const double *translationICP = ICP.GetTranslation();
325  const double *rotationDesign = Mtarget.GetRotationMatrix();
326  const double *rotationICP = ICP.GetRotationMatrix();
327 
328  double dx = translationDesign[0] - translationICP[0];
329  double dy = translationDesign[1] - translationICP[1];
330  double da = rotationDesign[1] - rotationICP[1];
331 
332  result.push_back(dx);
333  result.push_back(dy);
334  result.push_back(da);
335 
336  return result;
337 }
double dy
static TGeoHMatrix readTGeoHMatrix(std::string filename)
TGeoHMatrix getOverlapMatrixLikeICP(PndLmdOverlapInfo &info)
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
double dx
TGeoHMatrix PndLmdAlignQA::getMatrixSensorToSensor ( int  sensorOne,
int  sensorTwo 
)

Definition at line 177 of file PndLmdAlignQA.cxx.

References PndLmdGeometryHelper::getMatrixPndGlobalToSensor(), and helper.

Referenced by calculateOverlappingArea(), checkCombined(), getMisalignedOverlapFromGeoManager(), getMisalignedOverlapFromICP(), and getOverlapMatrixLikeICP().

177  {
178 
179  // the passive transformation that base-transforms a vector in Panda Global to the system of Sensor x
180  TGeoHMatrix toSen0 = helper->getMatrixPndGlobalToSensor(sensorOne);
181  TGeoHMatrix toSen1 = helper->getMatrixPndGlobalToSensor(sensorTwo);
182 
183  return toSen0.Inverse() * toSen1;
184 
185 }
const TGeoHMatrix getMatrixPndGlobalToSensor(const int sensorId)
PndLmdGeometryHelper * helper
Definition: PndLmdAlignQA.h:45
TGeoHMatrix PndLmdAlignQA::getMisalignedOverlapFromGeoManager ( PndLmdOverlapInfo info)

Definition at line 159 of file PndLmdAlignQA.cxx.

References getMatrixSensorToSensor(), getOverlapMatrixLikeICP(), PndLmdOverlapInfo::id1, and PndLmdOverlapInfo::id2.

Referenced by checkCombined().

159  {
160 
161  TGeoHMatrix sen0to1 = getMatrixSensorToSensor(info.id1, info.id2);
162  TGeoHMatrix sen0to1ICPcorr = getOverlapMatrixLikeICP(info);
163  // this order is important!
164  return sen0to1ICPcorr * sen0to1;
165 }
TGeoHMatrix getOverlapMatrixLikeICP(PndLmdOverlapInfo &info)
TGeoHMatrix getMatrixSensorToSensor(int sensorOne, int sensorTwo)
TGeoHMatrix PndLmdAlignQA::getMisalignedOverlapFromICP ( PndLmdOverlapInfo info,
std::string  ICPmatrix 
)

Definition at line 168 of file PndLmdAlignQA.cxx.

References getMatrixSensorToSensor(), PndLmdOverlapInfo::id1, PndLmdOverlapInfo::id2, manager, and PndLmdAlignManager::readTGeoHMatrix().

Referenced by checkCombined().

168  {
169 
170  TGeoHMatrix sen0to1 = getMatrixSensorToSensor(info.id1, info.id2);
171  TGeoHMatrix sen0to1ICPcorr = manager.readTGeoHMatrix(ICPmatrix);
172  // this order is important!
173  return sen0to1ICPcorr * sen0to1;
174 }
static TGeoHMatrix readTGeoHMatrix(std::string filename)
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
TGeoHMatrix getMatrixSensorToSensor(int sensorOne, int sensorTwo)
TGeoHMatrix PndLmdAlignQA::getOverlapMatrixLikeICP ( PndLmdOverlapInfo info)

Definition at line 187 of file PndLmdAlignQA.cxx.

References baseTransformation(), getMatrixSensorToSensor(), PndLmdOverlapInfo::id1, PndLmdOverlapInfo::id2, matricesMisaligned, PndLmdOverlapInfo::path1, and PndLmdOverlapInfo::path2.

Referenced by getMatrixDiffCM(), and getMisalignedOverlapFromGeoManager().

187  {
188 
189  if (!matricesMisaligned) {
190  cerr << "Error! No matrices in memory, please load them first.\n";
191  return TGeoHMatrix();
192  }
193 
194  //prepare misalignment matrices
195  TGeoHMatrix misalignmentToSensor0 = (*matricesMisaligned)[info.path1];
196  TGeoHMatrix misalignmentToSensor1 = (*matricesMisaligned)[info.path2];
197  // transform second misaignment Matrix to new base
198  TGeoHMatrix Sen0ToSen1 = getMatrixSensorToSensor(info.id1, info.id2);
199  TGeoHMatrix misSen1inSen0 = baseTransformation(misalignmentToSensor1, Sen0ToSen1);
200 
201  return misalignmentToSensor0.Inverse() * misSen1inSen0;
202 }
std::map< std::string, TGeoHMatrix > * matricesMisaligned
Definition: PndLmdAlignQA.h:51
TGeoHMatrix getMatrixSensorToSensor(int sensorOne, int sensorTwo)
TGeoHMatrix baseTransformation(TGeoHMatrix &input, TGeoHMatrix &toBaseMatrix)
PndLmdOverlapInfo & PndLmdAlignQA::getSmallOverlapInfo ( std::vector< PndLmdOverlapInfo > &  infos,
int  smallOverlap 
)

Definition at line 339 of file PndLmdAlignQA.cxx.

References PndLmdOverlapInfo::overlapID.

Referenced by checkCombined().

340  {
341  //PndLmdOverlapInfo result;
342  for (auto &info : infos) {
343  int smallOverlapHere = info.overlapID % 10;
344  if (smallOverlap == smallOverlapHere) {
345  return info;
346  }
347  }
348  return infos[0];
349 }
void PndLmdAlignQA::init ( )

Definition at line 44 of file PndLmdAlignQA.cxx.

References _inCentimeters, byPlane, PndLmdGeometryHelper::getInstance(), helper, infoAbsolute, infoMomentum, and infoRelative.

Referenced by PndLmdAlignQA().

44  {
45  infoMomentum = -1;
46  infoAbsolute = false;
47  infoRelative = false;
48  byPlane = false;
49  _inCentimeters = false;
51 }
double infoMomentum
Definition: PndLmdAlignQA.h:40
static PndLmdGeometryHelper & getInstance()
PndLmdGeometryHelper * helper
Definition: PndLmdAlignQA.h:45
int PndLmdAlignQA::noOfPairs ( int  overlapID)
private

Definition at line 274 of file PndLmdAlignQA.cxx.

References matrixInfo.

274  {
275  return matrixInfo[overlapID];
276 }
std::map< int, int > matrixInfo
Definition: PndLmdAlignQA.h:48
void PndLmdAlignQA::readMatrixInfo ( )

Definition at line 88 of file PndLmdAlignQA.cxx.

References _inCentimeters, filename, PndLmdAlignManager::findRegex(), LMDMatPath, manager, and matrixInfo.

Referenced by checkCyclicMatrices().

88  {
89 
90  string filename;
91  if (_inCentimeters) {
92  filename = LMDMatPath + "/info-cm.txt";
93  }
94  else {
95  filename = LMDMatPath + "/info-px.txt";
96  }
97 
98  std::ifstream info;
99  info.open(filename.c_str());
100  if (info.is_open()) {
101 
102  string line;
103  std::vector<string> values;
104  int overlapid, noPairs;
105  bool alignerComplete = true;
106 
107  while (std::getline(info, line)) {
108  if (line.find("aligner") != std::string::npos && alignerComplete) {
109  //cout << "found aligner line\n";
110  values = manager.findRegex(line, "aligner (\\d{1,4})");
111  if (values.size() > 1) {
112  overlapid = std::stoi(values[1]);
113  //cout << "aligner: " << overlapid << endl;
114  alignerComplete = false;
115  }
116  }
117  else if (line.find("no of pairs") != std::string::npos) {
118  //cout << "found pairs line\n";
119  values = manager.findRegex(line, "no of pairs. (\\d{1,6})");
120  if (values.size() > 1) {
121  noPairs = std::stoi(values[1]);
122  matrixInfo[overlapid] = noPairs;
123  alignerComplete = true;
124  }
125  }
126  }
127  }
128  else {
129  throw std::runtime_error("Unable to open file " + filename);
130  }
131 }
PndLmdAlignManager manager
Definition: PndLmdAlignQA.h:44
std::string LMDMatPath
Definition: PndLmdAlignQA.h:39
std::map< int, int > matrixInfo
Definition: PndLmdAlignQA.h:48
static std::vector< std::string > findRegex(std::string source, std::string regex)
const string filename
std::map< std::string, TGeoHMatrix > * PndLmdAlignQA::readRootMatrices ( TString filename)

Definition at line 133 of file PndLmdAlignQA.cxx.

Referenced by checkCombined().

133  {
134  cout << "reading matrices from file: " << filename << "\n";
135 
136  TFile *misalignmentMatrixRootfile = new TFile(filename, "READ");
137 
138  if (misalignmentMatrixRootfile->IsOpen()) {
139  std::map<std::string, TGeoHMatrix> *matrices;
140 
141  gDirectory->GetObject("PndLmdMisalignMatrices", matrices);
142  misalignmentMatrixRootfile->Close();
143  cout << matrices->size() << " matrices successfully read from file.\n";
144 
145  return matrices;
146  }
147  else {
148  cout << "file could not be read\n";
149  return NULL;
150  }
151 
152 }
const string filename
void PndLmdAlignQA::setAlignedGeometry ( bool  aligned)
inline

Definition at line 92 of file PndLmdAlignQA.h.

92 { alignOptionBool = aligned; }
bool alignOptionBool
Definition: PndLmdAlignQA.h:42
void PndLmdAlignQA::setBinaryMatPath ( const std::string &  path)
inline

Definition at line 90 of file PndLmdAlignQA.h.

90 { binaryMatPath = path; }
std::string binaryMatPath
Definition: PndLmdAlignQA.h:39
void PndLmdAlignQA::setInCentimeters ( bool  inCentimeters)
inline

Definition at line 88 of file PndLmdAlignQA.h.

Referenced by runLumiPixel2gAlignQA().

88 { this->_inCentimeters = inCentimeters; }
void PndLmdAlignQA::setInfoAbsolute ( bool  info)
inline

Definition at line 95 of file PndLmdAlignQA.h.

95 { infoAbsolute = info; }
void PndLmdAlignQA::setInfoMomentum ( double  info)
inline

Definition at line 96 of file PndLmdAlignQA.h.

96 { infoMomentum = info; }
double infoMomentum
Definition: PndLmdAlignQA.h:40
void PndLmdAlignQA::setInfoRelative ( bool  info)
inline

Definition at line 97 of file PndLmdAlignQA.h.

97 { infoRelative = info; }
void PndLmdAlignQA::setLmdMatPath ( const std::string &  path)
inline

Definition at line 89 of file PndLmdAlignQA.h.

Referenced by runLumiPixel2gAlignQA().

89 { LMDMatPath = path; }
std::string LMDMatPath
Definition: PndLmdAlignQA.h:39
void PndLmdAlignQA::setOutputPath ( std::string  path)
inline

Definition at line 94 of file PndLmdAlignQA.h.

94 { outputPath = path; }
std::string outputPath
Definition: PndLmdAlignQA.h:39
void PndLmdAlignQA::setPairsRequired ( int  number)
inline

Definition at line 93 of file PndLmdAlignQA.h.

Referenced by runLumiPixel2gAlignQA().

93 { pairsRequired = number; }
void PndLmdAlignQA::setPdfOutPath ( const std::string &  path)
inline

Definition at line 91 of file PndLmdAlignQA.h.

Referenced by runLumiPixel2gAlignQA().

91 { pdfOutPath = path; }
std::string pdfOutPath
Definition: PndLmdAlignQA.h:39

Member Data Documentation

bool PndLmdAlignQA::_inCentimeters
private
std::vector<std::string> PndLmdAlignQA::_inputFiles
private

Definition at line 38 of file PndLmdAlignQA.h.

bool PndLmdAlignQA::alignOptionBool
private

Definition at line 42 of file PndLmdAlignQA.h.

std::string PndLmdAlignQA::binaryMatPath
private

Definition at line 39 of file PndLmdAlignQA.h.

bool PndLmdAlignQA::byPlane
private

Definition at line 41 of file PndLmdAlignQA.h.

Referenced by init().

PndLmdGeometryHelper* PndLmdAlignQA::helper
private
bool PndLmdAlignQA::infoAbsolute
private

Definition at line 41 of file PndLmdAlignQA.h.

Referenced by init().

double PndLmdAlignQA::infoMomentum
private

Definition at line 40 of file PndLmdAlignQA.h.

Referenced by init().

bool PndLmdAlignQA::infoRelative
private

Definition at line 41 of file PndLmdAlignQA.h.

Referenced by init().

std::string PndLmdAlignQA::LMDMatPath
private

Definition at line 39 of file PndLmdAlignQA.h.

Referenced by checkCyclicMatrices(), checkForMatrixFiles(), and readMatrixInfo().

PndLmdAlignManager PndLmdAlignQA::manager
private
std::map<std::string, TGeoHMatrix>* PndLmdAlignQA::matricesMisaligned
private

Definition at line 51 of file PndLmdAlignQA.h.

Referenced by checkCombined(), and getOverlapMatrixLikeICP().

std::map<int, int> PndLmdAlignQA::matrixInfo
private

Definition at line 48 of file PndLmdAlignQA.h.

Referenced by noOfPairs(), and readMatrixInfo().

std::string PndLmdAlignQA::outputPath
private

Definition at line 39 of file PndLmdAlignQA.h.

int PndLmdAlignQA::pairsRequired
private

Definition at line 43 of file PndLmdAlignQA.h.

std::string PndLmdAlignQA::pdfOutPath
private

Definition at line 39 of file PndLmdAlignQA.h.

Referenced by checkCombined(), and checkCyclicMatrices().


The documentation for this class was generated from the following files: