Belle II Software  release-08-01-10
InverseRaytracer Class Reference

Utility for solving inverse ray-tracing problem. More...

#include <InverseRaytracer.h>

Inheritance diagram for InverseRaytracer:
Collaboration diagram for InverseRaytracer:

Classes

struct  CerenkovAngle
 Sine and cosine of Cerenkov angle. More...
 
struct  Solution
 Solution of inverse ray-tracing. More...
 

Public Types

enum  ErrorCodes {
  c_NoPhysicsSolution = -1 ,
  c_NoEquationSolution = -2 ,
  c_NoReflectionPoint = -3
}
 Error codes returned by solveDirect or solveReflected. More...
 
enum  EGeometry {
  c_Unified = 0 ,
  c_Segmented = 1
}
 Treatement of quartz geometry. More...
 
enum  EOptics {
  c_SemiLinear = 0 ,
  c_Exact = 1
}
 Treatement of spherical mirror optics. More...
 

Public Member Functions

 InverseRaytracer (int moduleID, double cosTotal)
 Class constructor. More...
 
void clear () const
 Clear the solutions to prepare for the new round.
 
int solveDirect (double xD, double zD, const TOPTrack::AssumedEmission &assumedEmission, const CerenkovAngle &cer, double step=0) const
 Solve inverse ray-tracing for direct photon. More...
 
int solveReflected (double xD, double zD, int Nxm, double xmMin, double xmMax, const TOPTrack::AssumedEmission &assumedEmission, const CerenkovAngle &cer, double step=0) const
 Solve inverse ray-tracing for reflected photon. More...
 
int solveForReflectionPoint (double xM, int Nxm, const TOPTrack::AssumedEmission &assumedEmission, const CerenkovAngle &cer) const
 Solve inverse ray-tracing for a given reflection point on the mirror. More...
 
double getCosTotal () const
 Returns cosine of total reflection angle. More...
 
std::vector< Solution > & getSolutions (unsigned i) const
 Returns the solutions of inverse ray-tracing. More...
 
bool getStatus () const
 Returns status. More...
 
bool getStatus (unsigned i) const
 Returns status. More...
 
bool isNymDifferent () const
 Checks if Nym differs between solutions front and back in at least one of the vectors. More...
 
PhotonState getReconstructedPhoton (const Solution &sol, double DFic=0) const
 Returns reconstructed photon at emission for a given solution of inverse raytracing. More...
 
int getModuleID () const
 Returns slot ID. More...
 
EGeometry getGeometry () const
 Returns quartz geometry treatement. More...
 
EOptics getOptics () const
 Returns treatement of spherical mirror optics. More...
 
const std::vector< BarSegment > & getBars () const
 Returns geometry data of bar segments. More...
 
const MirrorgetMirror () const
 Returns geometry data of spherical mirror. More...
 
const PrismgetPrism () const
 Returns geometry data of prism. More...
 
void setMirrorCenter (double xc, double yc)
 Sets the mirror center-of-curvature. More...
 

Static Public Member Functions

static void setMaxPropagationLen (double maxLen)
 Sets maximal allowed propagation length. More...
 

Protected Attributes

int m_moduleID = 0
 slot ID
 
EGeometry m_geometry = c_Unified
 quartz geometry
 
EOptics m_optics = c_SemiLinear
 spherical mirror optics
 
std::vector< BarSegmentm_bars
 geometry data of bar segments
 
Mirror m_mirror
 spherical mirror geometry data
 
Prism m_prism
 prism geometry data
 

Private Member Functions

bool solve (double dxdz, const CerenkovAngle &cer, const TOPTrack::TrackAngles &trk) const
 Solve inverse ray-tracing for unknown azimuthal Cerenkov angle. More...
 
bool findReflectionPoint (double xE, double zE, double xD, double zD, double xmMin, double xmMax, double &xM, double &zM, double &dxdz) const
 Finds reflection point on the mirror using semi-linear optics approximation. More...
 
double getDeltaXE (double x, double xe, double ze, double xd, double zd) const
 Returns the difference between input xe and the reflected position at given x. More...
 

Private Attributes

double m_cosTotal = 0
 cosine of total reflection angle
 
std::vector< Solutionm_solutions [2]
 storage for the two solutions
 
bool m_ok [2] = {false, false}
 status of solutions
 
ROOT::Math::XYZPoint m_emiPoint
 temporary storage of emission point
 
CerenkovAngle m_cer
 temporary storage of Cerenkov angle
 
TOPTrack::TrackAngles m_trk
 temporary storage of track polar and azimuthal angles
 

Static Private Attributes

static double s_maxLen = 10000
 maximal allowed propagation length
 

Detailed Description

Utility for solving inverse ray-tracing problem.

Definition at line 28 of file InverseRaytracer.h.

Member Enumeration Documentation

◆ EGeometry

enum EGeometry
inherited

Treatement of quartz geometry.

Enumerator
c_Unified 

single bar with average width and thickness

c_Segmented 

segmented bars

Definition at line 33 of file RaytracerBase.h.

◆ EOptics

enum EOptics
inherited

Treatement of spherical mirror optics.

Enumerator
c_SemiLinear 

semi-linear approximation

c_Exact 

exact optics

Definition at line 41 of file RaytracerBase.h.

◆ ErrorCodes

enum ErrorCodes

Error codes returned by solveDirect or solveReflected.

Enumerator
c_NoPhysicsSolution 

no physics solution

c_NoEquationSolution 

no solution of equation

c_NoReflectionPoint 

position on the mirror not found

Definition at line 35 of file InverseRaytracer.h.

35  {
36  c_NoPhysicsSolution = -1,
39  };
@ c_NoEquationSolution
no solution of equation
@ c_NoReflectionPoint
position on the mirror not found
@ c_NoPhysicsSolution
no physics solution

Constructor & Destructor Documentation

◆ InverseRaytracer()

InverseRaytracer ( int  moduleID,
double  cosTotal 
)
inline

Class constructor.

Parameters
moduleIDslot ID
cosTotalcosine of total reflection angle

Definition at line 133 of file InverseRaytracer.h.

Member Function Documentation

◆ findReflectionPoint()

bool findReflectionPoint ( double  xE,
double  zE,
double  xD,
double  zD,
double  xmMin,
double  xmMax,
double &  xM,
double &  zM,
double &  dxdz 
) const
private

Finds reflection point on the mirror using semi-linear optics approximation.

Parameters
xEunfolded emission position in x (unfolding w.r.t mirror)
zEemission position in z
xDunfolded detection position in x (unfolding w.r.t mirror)
zDdetection position in z
xmMinlower limit for search range in x
xmMaxupper limit for search range in x
xMreflection position in x [out]
zMreflection position in z [out]
dxdzphoton slope in x-z (dx/dz) before mirror [out]
Returns
true on success

Definition at line 313 of file InverseRaytracer.cc.

316  {
317  double xe = (xE - m_mirror.xc) / m_mirror.R;
318  double ze = (zE - m_mirror.zc) / m_mirror.R;
319  double xd = (xD - m_mirror.xc) / m_mirror.R;
320  double zd = (zD - m_mirror.zc) / m_mirror.R;
321 
322  double x1 = (xmMin - m_mirror.xc) / m_mirror.R;
323  double y1 = getDeltaXE(x1, xe, ze, xd, zd);
324 
325  double x2 = (xmMax - m_mirror.xc) / m_mirror.R;
326  double y2 = getDeltaXE(x2, xe, ze, xd, zd);
327 
328  if (y1 * y2 > 0) return false; // no (single) solution
329 
330  for (int i = 0; i < 20; i++) {
331  double x = (x1 + x2) / 2;
332  double y = getDeltaXE(x, xe, ze, xd, zd);
333  if (y * y1 < 0) {
334  x2 = x;
335  } else {
336  x1 = x;
337  y1 = y;
338  }
339  }
340  double x = (x1 + x2) / 2;
341  double z = sqrt(1 - x * x);
342  xM = x * m_mirror.R + m_mirror.xc;
343  zM = z * m_mirror.R + m_mirror.zc;
344 
345  double kx = (x - xd);
346  double kz = (z - zd);
347  double s = 2 * (kx * x + kz * z);
348  kx -= s * x;
349  kz -= s * z;
350  dxdz = kx / kz;
351 
352  return true;
353  }
double getDeltaXE(double x, double xe, double ze, double xd, double zd) const
Returns the difference between input xe and the reflected position at given x.
Mirror m_mirror
spherical mirror geometry data
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
double xc
center of curvature in x
Definition: RaytracerBase.h:81
double zc
center of curvature in z
Definition: RaytracerBase.h:83

◆ getBars()

const std::vector<BarSegment>& getBars ( ) const
inlineinherited

Returns geometry data of bar segments.

Returns
geometry data of bar segments

Definition at line 161 of file RaytracerBase.h.

◆ getCosTotal()

double getCosTotal ( ) const
inline

Returns cosine of total reflection angle.

Returns
cosine of total reflection angle at mean photon energy for beta = 1

Definition at line 191 of file InverseRaytracer.h.

◆ getDeltaXE()

double getDeltaXE ( double  x,
double  xe,
double  ze,
double  xd,
double  zd 
) const
private

Returns the difference between input xe and the reflected position at given x.

This function is used to find the reflection point on the mirror. All arguments must be given in the mirror frame and in units of mirror radius.

Parameters
xposition in x on the mirror
xeunfolded emission position in x
zeemission position in z
xdunfolded detection position in x
zddetection position in z
Returns
the difference in units of mirror radius.

Definition at line 300 of file InverseRaytracer.cc.

◆ getGeometry()

EGeometry getGeometry ( ) const
inlineinherited

Returns quartz geometry treatement.

Returns
quartz geometry treatement

Definition at line 149 of file RaytracerBase.h.

◆ getMirror()

const Mirror& getMirror ( ) const
inlineinherited

Returns geometry data of spherical mirror.

Returns
geometry data of spherical mirror

Definition at line 167 of file RaytracerBase.h.

◆ getModuleID()

int getModuleID ( ) const
inlineinherited

Returns slot ID.

Returns
slot ID

Definition at line 143 of file RaytracerBase.h.

◆ getOptics()

EOptics getOptics ( ) const
inlineinherited

Returns treatement of spherical mirror optics.

Returns
spherical mirror optics

Definition at line 155 of file RaytracerBase.h.

◆ getPrism()

const Prism& getPrism ( ) const
inlineinherited

Returns geometry data of prism.

Returns
geometry data of prism

Definition at line 173 of file RaytracerBase.h.

◆ getReconstructedPhoton()

PhotonState getReconstructedPhoton ( const Solution sol,
double  DFic = 0 
) const

Returns reconstructed photon at emission for a given solution of inverse raytracing.

Parameters
solsolution of inverse raytracing
DFicadditional rotation angle around track (delta Cerenkov azimuthal angle)
Returns
reconstructed photon (note: usually meaningless if status is false!)

Definition at line 356 of file InverseRaytracer.cc.

◆ getSolutions()

std::vector<Solution>& getSolutions ( unsigned  i) const
inline

Returns the solutions of inverse ray-tracing.

Parameters
ito select first (0) or second (1) solutions
Returns
a vector of solutions

Definition at line 198 of file InverseRaytracer.h.

◆ getStatus() [1/2]

bool getStatus ( void  ) const
inline

Returns status.

Returns
status

Definition at line 204 of file InverseRaytracer.h.

◆ getStatus() [2/2]

bool getStatus ( unsigned  i) const
inline

Returns status.

Parameters
ito select status of first (0) or second (1) solutions
Returns
status

Definition at line 211 of file InverseRaytracer.h.

◆ isNymDifferent()

bool isNymDifferent ( ) const

Checks if Nym differs between solutions front and back in at least one of the vectors.

Returns
true if differs

Definition at line 259 of file InverseRaytracer.cc.

◆ setMaxPropagationLen()

static void setMaxPropagationLen ( double  maxLen)
inlinestatic

Sets maximal allowed propagation length.

Parameters
maxLenmaximal allowed propagation length

Definition at line 141 of file InverseRaytracer.h.

◆ setMirrorCenter()

void setMirrorCenter ( double  xc,
double  yc 
)
inherited

Sets the mirror center-of-curvature.

Parameters
xccenter of curvature in x
yccenter of curvature in y

Definition at line 100 of file RaytracerBase.cc.

◆ solve()

bool solve ( double  dxdz,
const CerenkovAngle cer,
const TOPTrack::TrackAngles trk 
) const
private

Solve inverse ray-tracing for unknown azimuthal Cerenkov angle.

Solutions are always two or none (very rare case, phic=0: the two are the same)

Parameters
dxdzphoton slope (dx/dz) in x-z projection
cersine and cosine of Cerenkov angle
trksine and cosine of track polar and azimuthal angles at photon emission
Returns
true if solutions exits

Definition at line 268 of file InverseRaytracer.cc.

◆ solveDirect()

int solveDirect ( double  xD,
double  zD,
const TOPTrack::AssumedEmission assumedEmission,
const CerenkovAngle cer,
double  step = 0 
) const

Solve inverse ray-tracing for direct photon.

Parameters
xDunfolded position in x of photon at detection plane
zDposition of detection plane
assumedEmissionphoton emission position and track angles
cersine and cosine of Cerenkov angle
stepstep for numerical derivative calculation
Returns
index of solution or ErrorCode on fail

Definition at line 67 of file InverseRaytracer.cc.

◆ solveForReflectionPoint()

int solveForReflectionPoint ( double  xM,
int  Nxm,
const TOPTrack::AssumedEmission assumedEmission,
const CerenkovAngle cer 
) const

Solve inverse ray-tracing for a given reflection point on the mirror.

Parameters
xMposition in x of the reflection point on the mirror.
Nxmsigned number of reflections before mirror
assumedEmissionphoton emission position and track angles
cersine and cosine of Cerenkov angle
Returns
index of solution or ErrorCode on fail

Definition at line 192 of file InverseRaytracer.cc.

◆ solveReflected()

int solveReflected ( double  xD,
double  zD,
int  Nxm,
double  xmMin,
double  xmMax,
const TOPTrack::AssumedEmission assumedEmission,
const CerenkovAngle cer,
double  step = 0 
) const

Solve inverse ray-tracing for reflected photon.

Parameters
xDunfolded position in x of photon at detection plane (unfolding w.r.t mirror)
zDposition of detection plane
Nxmsigned number of reflections before mirror
xmMinlower limit for the reflection point search range
xmMaxupper limit for the reflection point search range
assumedEmissionphoton emission position and track angles
cersine and cosine of Cerenkov angle
stepstep for numerical derivative calculation
Returns
index of solution or ErrorCode on fail

Definition at line 108 of file InverseRaytracer.cc.


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