Project
Loading...
Searching...
No Matches
o2::mch::Aligner Class Reference

#include <Aligner.h>

Inherits TObject.

Public Types

enum  {
  fgNSt = 5 , fgNCh = 10 , fgNTrkMod = 16 , fgNHalfCh = 20 ,
  fgNDetHalfChMax = 13 , fgNDetElem = 156 , fNLocal = 4 , fgNParCh = 4 ,
  fNGlobal = fgNParCh * fgNDetElem
}
 
enum  ParameterMask {
  ParX = 1 << 0 , ParY = 1 << 1 , ParZ = 1 << 2 , ParTZ = 1 << 3 ,
  ParAllTranslations = ParX | ParY | ParZ , ParAllRotations = ParTZ , ParAll = ParAllTranslations | ParAllRotations
}
 global parameter bit set, used for masks More...
 
enum  SidesMask {
  SideTop = 1 << 0 , SideLeft = 1 << 1 , SideBottom = 1 << 2 , SideRight = 1 << 3 ,
  SideTopLeft = SideTop | SideLeft , SideTopRight = SideTop | SideRight , SideBottomLeft = SideBottom | SideLeft , SideBottomRight = SideBottom | SideRight ,
  AllSides = SideTop | SideBottom | SideLeft | SideRight
}
 detector sides bit set, used for selecting sides in constrains More...
 

Public Member Functions

 Aligner ()
 
 ~Aligner ()
 
void init (TString DataRecFName="millerecords.root", TString ConsRecFName="milleconstraints.root")
 
void terminate ()
 
void ProcessTrack (Track &track, const o2::mch::geo::TransformationCreator &transformation, Bool_t doAlignment, Double_t weight=1)
 
void SetRunNumber (int id)
 run number
 
void SetBFieldOn (bool value)
 Set flag for Magnetic field On/Off.
 
void SetDoEvaluation (bool value)
 set to true to do refit evaluation
 
void SetRefitStraightTracks (bool value)
 set to true to refit tracks
 
void SetAllowedVariation (int iPar, double value)
 
void SetSigmaXY (double sigmaX, double sigmaY)
 
void FixAll (unsigned int parameterMask=ParAll)
 
void FixChamber (int iCh, unsigned int parameterMask=ParAll)
 
void FixDetElem (int iDetElemId, unsigned int parameterMask=ParAll)
 
void FixHalfSpectrometer (const bool *bChOnOff, unsigned int sidesMask=AllSides, unsigned int parameterMask=ParAll)
 
void FixParameter (int iPar)
 
void FixParameter (int iDetElem, int iPar)
 
void ReleaseChamber (int iCh, unsigned int parameterMask=ParAll)
 
void ReleaseDetElem (int iDetElemId, unsigned int parameterMask=ParAll)
 
void ReleaseParameter (int iPar)
 
void ReleaseParameter (int iDetElem, int iPar)
 
void GroupChamber (int iCh, unsigned int parameterMask=ParAll)
 
void GroupHalfChamber (int iCh, int iHalf, unsigned int parameterMask=ParAll)
 
void GroupDetElems (int detElemMin, int detElemMax, unsigned int parameterMask=ParAll)
 
void GroupDetElems (const int *detElemList, int nDetElem, unsigned int parameterMask=ParAll)
 
void SetChamberNonLinear (int iCh, unsigned int parameterMask)
 
void SetDetElemNonLinear (int iSt, unsigned int parameterMask)
 
void SetParameterNonLinear (int iPar)
 
void SetParameterNonLinear (int iDetElem, int iPar)
 
void AddConstraints (const bool *bChOnOff, unsigned int parameterMask)
 
void AddConstraints (const bool *bChOnOff, const bool *lVarXYT, unsigned int sidesMask=AllSides)
 
void InitGlobalParameters (double *par)
 initialize global parameters to a give set of values
 
void GlobalFit (std::vector< double > &params, std::vector< double > &errors, std::vector< double > &pulls)
 perform global fit
 
void PrintGlobalParameters (void) const
 print global parameters
 
double GetParError (int iPar) const
 get error on a given parameter
 
o2::fwdalign::MillePedeRecordGetRecord ()
 
void ReAlign (std::vector< o2::detectors::AlignParam > &params, std::vector< double > &misAlignments)
 
void SetAlignmentResolution (const TClonesArray *misAlignArray, int chId, double chResX, double chResY, double deResX, double deResY)
 
TTree * GetResTree ()
 
void SetReadOnly ()
 
void DisableRecordWriter ()
 

Static Public Attributes

static const int fgNDetElemCh [fgNCh] = {4, 4, 4, 4, 18, 18, 26, 26, 26, 26}
 Number of detection elements per chamber.
 
static const int fgSNDetElemCh [fgNCh+1] = {0, 4, 8, 12, 16, 34, 52, 78, 104, 130, 156}
 Sum of detection elements up to this chamber.
 
static const int fgNDetElemHalfCh [fgNHalfCh] = {2, 2, 2, 2, 2, 2, 2, 2, 9, 9, 9, 9, 13, 13, 13, 13, 13, 13, 13, 13}
 Number of detection element per tracking module.
 
static const int fgDetElemHalfCh [fgNHalfCh][fgNDetHalfChMax]
 list of detection elements per tracking module
 

Detailed Description

Definition at line 104 of file Aligner.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
fgNSt 

Number tracking stations.

fgNCh 

Number tracking chambers.

fgNTrkMod 

Number of tracking modules.

fgNHalfCh 

Number of half chambers.

fgNDetHalfChMax 

max number of detector elements per half chamber

fgNDetElem 

Total number of detection elements (4*2 + 4*2 + 18*2 + 26*2 + 26*2)

fNLocal 

Number of local parameters.

fgNParCh 

Number of degrees of freedom per chamber.

fNGlobal 

Number of global parameters.

Definition at line 119 of file Aligner.h.

◆ ParameterMask

global parameter bit set, used for masks

Enumerator
ParX 
ParY 
ParZ 
ParTZ 
ParAllTranslations 
ParAllRotations 
ParAll 

Definition at line 162 of file Aligner.h.

◆ SidesMask

detector sides bit set, used for selecting sides in constrains

Enumerator
SideTop 
SideLeft 
SideBottom 
SideRight 
SideTopLeft 
SideTopRight 
SideBottomLeft 
SideBottomRight 
AllSides 

Definition at line 175 of file Aligner.h.

Constructor & Destructor Documentation

◆ Aligner()

o2::mch::Aligner::Aligner ( )

constructor

Definition at line 126 of file Aligner.cxx.

◆ ~Aligner()

o2::mch::Aligner::~Aligner ( )

Definition at line 187 of file Aligner.cxx.

Member Function Documentation

◆ AddConstraints() [1/2]

void o2::mch::Aligner::AddConstraints ( const bool *  bChOnOff,
const bool *  lVarXYT,
unsigned int  sidesMask = AllSides 
)

Add constraint equations for selected chambers, degrees of freedom and detector half

Definition at line 954 of file Aligner.cxx.

◆ AddConstraints() [2/2]

void o2::mch::Aligner::AddConstraints ( const bool *  bChOnOff,
unsigned int  parameterMask 
)

Add constraint equations for selected chambers and degrees of freedom

Definition at line 909 of file Aligner.cxx.

◆ DisableRecordWriter()

void o2::mch::Aligner::DisableRecordWriter ( )
inline

Definition at line 319 of file Aligner.h.

◆ FixAll()

void o2::mch::Aligner::FixAll ( unsigned int  parameterMask = ParAll)

fix parameters matching mask, for all chambers

Definition at line 554 of file Aligner.cxx.

◆ FixChamber()

void o2::mch::Aligner::FixChamber ( int  iCh,
unsigned int  parameterMask = ParAll 
)

fix parameters matching mask, for all detector elements in a given chamber, counting from 1

Definition at line 577 of file Aligner.cxx.

◆ FixDetElem()

void o2::mch::Aligner::FixDetElem ( int  iDetElemId,
unsigned int  parameterMask = ParAll 
)

fix parameters matching mask, for a given detector element, counting from 0

Definition at line 609 of file Aligner.cxx.

◆ FixHalfSpectrometer()

void o2::mch::Aligner::FixHalfSpectrometer ( const bool *  bChOnOff,
unsigned int  sidesMask = AllSides,
unsigned int  parameterMask = ParAll 
)

Fix parameters matching mask for all detectors in selected chambers and selected sides of the spectrometer

Definition at line 628 of file Aligner.cxx.

◆ FixParameter() [1/2]

void o2::mch::Aligner::FixParameter ( int  iDetElem,
int  iPar 
)
inline

Definition at line 230 of file Aligner.h.

◆ FixParameter() [2/2]

void o2::mch::Aligner::FixParameter ( int  iPar)

fix a given parameter, counting from 0

Definition at line 698 of file Aligner.cxx.

◆ GetParError()

double o2::mch::Aligner::GetParError ( int  iPar) const

get error on a given parameter

Definition at line 1316 of file Aligner.cxx.

◆ GetRecord()

o2::fwdalign::MillePedeRecord & o2::mch::Aligner::GetRecord ( )
inline

Definition at line 303 of file Aligner.h.

◆ GetResTree()

TTree * o2::mch::Aligner::GetResTree ( )
inline

Definition at line 309 of file Aligner.h.

◆ GlobalFit()

void o2::mch::Aligner::GlobalFit ( std::vector< double > &  params,
std::vector< double > &  errors,
std::vector< double > &  pulls 
)

perform global fit

Call global fit; Global parameters are stored in parameters

Definition at line 1298 of file Aligner.cxx.

◆ GroupChamber()

void o2::mch::Aligner::GroupChamber ( int  iCh,
unsigned int  parameterMask = ParAll 
)

group parameters matching mask for all detector elements in a given chamber, counting from 1

Definition at line 773 of file Aligner.cxx.

◆ GroupDetElems() [1/2]

void o2::mch::Aligner::GroupDetElems ( const int detElemList,
int  nDetElem,
unsigned int  parameterMask = ParAll 
)

group parameters matching mask for all detector elements in list

Definition at line 823 of file Aligner.cxx.

◆ GroupDetElems() [2/2]

void o2::mch::Aligner::GroupDetElems ( int  detElemMin,
int  detElemMax,
unsigned int  parameterMask = ParAll 
)

group parameters matching mask for all detector elements between min and max

Definition at line 802 of file Aligner.cxx.

◆ GroupHalfChamber()

void o2::mch::Aligner::GroupHalfChamber ( int  iCh,
int  iHalf,
unsigned int  parameterMask = ParAll 
)

group parameters matching mask for all detector elements in a given tracking module (= half chamber), counting from 0

Definition at line 786 of file Aligner.cxx.

◆ init()

void o2::mch::Aligner::init ( TString  DataRecFName = "millerecords.root",
TString  ConsRecFName = "milleconstraints.root" 
)

initialize

initialize millipede must be called after necessary detectors have been fixed, but before constrains are added and before global parameters initial value are set

Definition at line 195 of file Aligner.cxx.

◆ InitGlobalParameters()

void o2::mch::Aligner::InitGlobalParameters ( double *  par)

initialize global parameters to a give set of values

Initialize global parameters with par array

Definition at line 1256 of file Aligner.cxx.

◆ PrintGlobalParameters()

void o2::mch::Aligner::PrintGlobalParameters ( void  ) const

print global parameters

Definition at line 1310 of file Aligner.cxx.

◆ ProcessTrack()

void o2::mch::Aligner::ProcessTrack ( Track track,
const o2::mch::geo::TransformationCreator transformation,
Bool_t  doAlignment,
Double_t  weight = 1 
)

process track for alignment minimization

Definition at line 376 of file Aligner.cxx.

◆ ReAlign()

void o2::mch::Aligner::ReAlign ( std::vector< o2::detectors::AlignParam > &  params,
std::vector< double > &  misAlignments 
)

Returns a new AliMUONGeometryTransformer with the found misalignments applied.

Definition at line 1322 of file Aligner.cxx.

◆ ReleaseChamber()

void o2::mch::Aligner::ReleaseChamber ( int  iCh,
unsigned int  parameterMask = ParAll 
)

release parameters matching mask, for all detector elements in a given chamber, counting from 1

Definition at line 710 of file Aligner.cxx.

◆ ReleaseDetElem()

void o2::mch::Aligner::ReleaseDetElem ( int  iDetElemId,
unsigned int  parameterMask = ParAll 
)

release parameters matching mask, for a given detector element, counting from 0

Definition at line 742 of file Aligner.cxx.

◆ ReleaseParameter() [1/2]

void o2::mch::Aligner::ReleaseParameter ( int  iDetElem,
int  iPar 
)
inline

Definition at line 246 of file Aligner.h.

◆ ReleaseParameter() [2/2]

void o2::mch::Aligner::ReleaseParameter ( int  iPar)

release a given parameter, counting from 0

Definition at line 761 of file Aligner.cxx.

◆ SetAlignmentResolution()

void o2::mch::Aligner::SetAlignmentResolution ( const TClonesArray *  misAlignArray,
int  chId,
double  chResX,
double  chResY,
double  deResX,
double  deResY 
)

Set alignment resolution to misalign objects to be stored in CDB if rChId is > 0 set parameters for this chamber only, counting from 1

Definition at line 1404 of file Aligner.cxx.

◆ SetAllowedVariation()

void o2::mch::Aligner::SetAllowedVariation ( int  iPar,
double  value 
)

"Encouraged" variation for degrees of freedom

Definition at line 1267 of file Aligner.cxx.

◆ SetBFieldOn()

void o2::mch::Aligner::SetBFieldOn ( bool  value)
inline

Set flag for Magnetic field On/Off.

Definition at line 199 of file Aligner.h.

◆ SetChamberNonLinear()

void o2::mch::Aligner::SetChamberNonLinear ( int  iCh,
unsigned int  parameterMask 
)

Set parameters matching mask as non linear, for all detector elements in a given chamber, counting from 1

Definition at line 855 of file Aligner.cxx.

◆ SetDetElemNonLinear()

void o2::mch::Aligner::SetDetElemNonLinear ( int  iSt,
unsigned int  parameterMask 
)

Set parameters matching mask as non linear, for a given detector element, counting from 0

Definition at line 878 of file Aligner.cxx.

◆ SetDoEvaluation()

void o2::mch::Aligner::SetDoEvaluation ( bool  value)
inline

set to true to do refit evaluation

Definition at line 205 of file Aligner.h.

◆ SetParameterNonLinear() [1/2]

void o2::mch::Aligner::SetParameterNonLinear ( int  iDetElem,
int  iPar 
)
inline

Definition at line 275 of file Aligner.h.

◆ SetParameterNonLinear() [2/2]

void o2::mch::Aligner::SetParameterNonLinear ( int  iPar)

Set nonlinear flag for parameter iPar

Definition at line 897 of file Aligner.cxx.

◆ SetReadOnly()

void o2::mch::Aligner::SetReadOnly ( )
inline

Definition at line 314 of file Aligner.h.

◆ SetRefitStraightTracks()

void o2::mch::Aligner::SetRefitStraightTracks ( bool  value)
inline

set to true to refit tracks

Definition at line 211 of file Aligner.h.

◆ SetRunNumber()

void o2::mch::Aligner::SetRunNumber ( int  id)
inline

run number

Definition at line 193 of file Aligner.h.

◆ SetSigmaXY()

void o2::mch::Aligner::SetSigmaXY ( double  sigmaX,
double  sigmaY 
)

Set expected measurement resolution

Definition at line 1284 of file Aligner.cxx.

◆ terminate()

void o2::mch::Aligner::terminate ( )

Definition at line 359 of file Aligner.cxx.

Member Data Documentation

◆ fgDetElemHalfCh

const int o2::mch::Aligner::fgDetElemHalfCh
static

list of detection elements per tracking module

Definition at line 159 of file Aligner.h.

◆ fgNDetElemCh

const int o2::mch::Aligner::fgNDetElemCh = {4, 4, 4, 4, 18, 18, 26, 26, 26, 26}
static

Number of detection elements per chamber.

Definition at line 150 of file Aligner.h.

◆ fgNDetElemHalfCh

const int o2::mch::Aligner::fgNDetElemHalfCh = {2, 2, 2, 2, 2, 2, 2, 2, 9, 9, 9, 9, 13, 13, 13, 13, 13, 13, 13, 13}
static

Number of detection element per tracking module.

Definition at line 156 of file Aligner.h.

◆ fgSNDetElemCh

const int o2::mch::Aligner::fgSNDetElemCh = {0, 4, 8, 12, 16, 34, 52, 78, 104, 130, 156}
static

Sum of detection elements up to this chamber.

Definition at line 153 of file Aligner.h.


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