SND@LHC Software
Loading...
Searching...
No Matches
genfit::TGeoMaterialInterface Class Reference

AbsMaterialInterface implementation for use with ROOT's TGeoManager. More...

#include <TGeoMaterialInterface.h>

Inheritance diagram for genfit::TGeoMaterialInterface:
Collaboration diagram for genfit::TGeoMaterialInterface:

Public Member Functions

 TGeoMaterialInterface ()
 
virtual ~TGeoMaterialInterface ()
 
bool initTrack (double posX, double posY, double posZ, double dirX, double dirY, double dirZ)
 Initialize the navigator at given position and with given direction. Returns true if the volume changed.
 
void getMaterialParameters (double &density, double &Z, double &A, double &radiationLength, double &mEE)
 Get material parameters in current material.
 
void getMaterialParameters (MaterialProperties &parameters)
 
double findNextBoundary (const RKTrackRep *rep, const M1x7 &state7, double sMax, bool varField=true)
 Make a step (following the curvature) until step length sMax or the next boundary is reached. After making a step to a boundary, the position has to be beyond the boundary, i.e. the current material has to be that beyond the boundary. The actual step made is returned.
 
double findNextBoundaryAndStepStraight (double sMax)
 
 ClassDef (TGeoMaterialInterface, 1)
 
 TGeoMaterialInterface ()
 
virtual ~TGeoMaterialInterface ()
 
bool initTrack (double posX, double posY, double posZ, double dirX, double dirY, double dirZ)
 Initialize the navigator at given position and with given direction. Returns true if the volume changed.
 
void getMaterialParameters (double &density, double &Z, double &A, double &radiationLength, double &mEE)
 Get material parameters in current material.
 
void getMaterialParameters (MaterialProperties &parameters)
 
double findNextBoundary (const RKTrackRep *rep, const M1x7 &state7, double sMax, bool varField=true)
 Make a step (following the curvature) until step length sMax or the next boundary is reached. After making a step to a boundary, the position has to be beyond the boundary, i.e. the current material has to be that beyond the boundary. The actual step made is returned.
 
double findNextBoundaryAndStepStraight (double sMax)
 
 ClassDef (TGeoMaterialInterface, 1)
 
- Public Member Functions inherited from genfit::AbsMaterialInterface
 AbsMaterialInterface ()
 
virtual ~AbsMaterialInterface ()
 
 AbsMaterialInterface ()
 
virtual ~AbsMaterialInterface ()
 

Detailed Description

AbsMaterialInterface implementation for use with ROOT's TGeoManager.

Definition at line 35 of file TGeoMaterialInterface.h.

Constructor & Destructor Documentation

◆ TGeoMaterialInterface() [1/2]

genfit::TGeoMaterialInterface::TGeoMaterialInterface ( )
inline

Definition at line 39 of file TGeoMaterialInterface.h.

39{};

◆ ~TGeoMaterialInterface() [1/2]

virtual genfit::TGeoMaterialInterface::~TGeoMaterialInterface ( )
inlinevirtual

Definition at line 40 of file TGeoMaterialInterface.h.

40{;};

◆ TGeoMaterialInterface() [2/2]

genfit::TGeoMaterialInterface::TGeoMaterialInterface ( )
inline

Definition at line 39 of file TGeoMaterialInterface.h.

39{};

◆ ~TGeoMaterialInterface() [2/2]

virtual genfit::TGeoMaterialInterface::~TGeoMaterialInterface ( )
inlinevirtual

Definition at line 40 of file TGeoMaterialInterface.h.

40{;};

Member Function Documentation

◆ ClassDef() [1/2]

genfit::TGeoMaterialInterface::ClassDef ( TGeoMaterialInterface  ,
 
)

◆ ClassDef() [2/2]

genfit::TGeoMaterialInterface::ClassDef ( TGeoMaterialInterface  ,
 
)

◆ findNextBoundary() [1/2]

double genfit::TGeoMaterialInterface::findNextBoundary ( const RKTrackRep rep,
const M1x7 state7,
double  sMax,
bool  varField = true 
)
virtual

Make a step (following the curvature) until step length sMax or the next boundary is reached. After making a step to a boundary, the position has to be beyond the boundary, i.e. the current material has to be that beyond the boundary. The actual step made is returned.

Implements genfit::AbsMaterialInterface.

Definition at line 88 of file TGeoMaterialInterface.cc.

92{
93 const double delta(1.E-2); // cm, distance limit beneath which straight-line steps are taken.
94 const double epsilon(1.E-1); // cm, allowed upper bound on arch
95 // deviation from straight line
96
97 M1x3 SA;
98 M1x7 state7, oldState7;
99 memcpy(oldState7, stateOrig, sizeof(state7));
100
101 int stepSign(sMax < 0 ? -1 : 1);
102
103 double s = 0; // trajectory length to boundary
104
105 const unsigned maxIt = 300;
106 unsigned it = 0;
107
108 // Initialize the geometry to the current location (set by caller).
109 gGeoManager->FindNextBoundary(fabs(sMax) - s);
110 double safety = gGeoManager->GetSafeDistance();
111 double slDist = gGeoManager->GetStep();
112 double step = slDist;
113
114 while (1) {
115 if (++it > maxIt){
116 Exception exc("TGeoMaterialInterface::findNextBoundary ==> maximum number of iterations exceeded",__LINE__,__FILE__);
117 exc.setFatal();
118 throw exc;
119 }
120
121 // No boundary in sight?
122 if (s + safety > fabs(sMax)) {
123 if (debug)
124 std::cout << " next boundary is further away than sMax \n";
125 return stepSign*(s + safety); //sMax;
126 }
127
128 // Are we at the boundary?
129 if (slDist < delta) {
130 if (debug)
131 std::cout << " very close to the boundary -> return"
132 << " stepSign*(s + slDist) = "
133 << stepSign << "*(" << s + slDist << ")\n";
134 return stepSign*(s + slDist);
135 }
136
137 // We have to find whether there's any boundary on our path.
138
139 // Follow curved arch, then see if we may have missed a boundary.
140 // Always propagate complete way from original start to avoid
141 // inconsistent extrapolations.
142 memcpy(state7, stateOrig, sizeof(state7));
143 rep->RKPropagate(state7, NULL, SA, stepSign*(s + step), varField);
144
145 // Straight line distance² between extrapolation finish and
146 // the end of the previously determined safe segment.
147 double dist2 = (pow(state7[0] - oldState7[0], 2)
148 + pow(state7[1] - oldState7[1], 2)
149 + pow(state7[2] - oldState7[2], 2));
150 // Maximal lateral deviation².
151 double maxDeviation2 = 0.25*(step*step - dist2);
152
153 if (step > safety
154 && maxDeviation2 > epsilon*epsilon) {
155 // Need to take a shorter step to reliably estimate material,
156 // but only if we didn't move by safety. In order to avoid
157 // issues with roundoff we check 'step > safety' instead of
158 // 'step != safety'. If we moved by safety, there couldn't have
159 // been a boundary that we missed along the path, but we could
160 // be on the boundary.
161
162 // Take a shorter step, but never shorter than safety.
163 step = std::max(step / 2, safety);
164 } else {
165 gGeoManager->PushPoint();
166 bool volChanged = initTrack(state7[0], state7[1], state7[2],
167 stepSign*state7[3], stepSign*state7[4],
168 stepSign*state7[5]);
169
170 if (volChanged) {
171 // Move back to start.
172 gGeoManager->PopPoint();
173
174 // Extrapolation may not take the exact step length we asked
175 // for, so it can happen that a requested step < safety takes
176 // us across the boundary. This is then the best estimate we
177 // can get of the distance to the boundary with the stepper.
178 if (step <= safety)
179 return stepSign*(s + step);
180
181 // Volume changed during the extrapolation. Take a shorter
182 // step, but never shorter than safety.
183 step = std::max(step / 2, safety);
184 } else {
185 // we're in the new place, the step was safe, advance
186 s += step;
187
188 memcpy(oldState7, state7, sizeof(state7));
189 gGeoManager->PopDummy(); // Pop stack, but stay in place.
190
191 gGeoManager->FindNextBoundary(fabs(sMax) - s);
192 safety = gGeoManager->GetSafeDistance();
193 step = slDist = gGeoManager->GetStep();
194 }
195 }
196 }
197}
bool initTrack(double posX, double posY, double posZ, double dirX, double dirY, double dirZ)
Initialize the navigator at given position and with given direction. Returns true if the volume chang...
double M1x7[1 *7]
Definition RKTools.h:36
double M1x3[1 *3]
Definition RKTools.h:33

◆ findNextBoundary() [2/2]

double genfit::TGeoMaterialInterface::findNextBoundary ( const RKTrackRep rep,
const M1x7 state7,
double  sMax,
bool  varField = true 
)
virtual

Make a step (following the curvature) until step length sMax or the next boundary is reached. After making a step to a boundary, the position has to be beyond the boundary, i.e. the current material has to be that beyond the boundary. The actual step made is returned.

Implements genfit::AbsMaterialInterface.

◆ findNextBoundaryAndStepStraight() [1/2]

double genfit::TGeoMaterialInterface::findNextBoundaryAndStepStraight ( double  sMax)
virtual

Implements genfit::AbsMaterialInterface.

Definition at line 201 of file TGeoMaterialInterface.cc.

201 {
202
203 gGeoManager->FindNextBoundaryAndStep(sMax);
204 return gGeoManager->GetStep();
205
206}

◆ findNextBoundaryAndStepStraight() [2/2]

double genfit::TGeoMaterialInterface::findNextBoundaryAndStepStraight ( double  sMax)
virtual

◆ getMaterialParameters() [1/4]

void genfit::TGeoMaterialInterface::getMaterialParameters ( double &  density,
double &  Z,
double &  A,
double &  radiationLength,
double &  mEE 
)
virtual

Get material parameters in current material.

Implements genfit::AbsMaterialInterface.

Definition at line 56 of file TGeoMaterialInterface.cc.

60 {
61
62 TGeoMaterial* mat = gGeoManager->GetCurrentVolume()->GetMedium()->GetMaterial();
63
64 density = mat->GetDensity();
65 Z = mat->GetZ();
66 A = mat->GetA();
67 radiationLength = mat->GetRadLen();
68 mEE = MeanExcEnergy_get(mat);
69
70}
double MeanExcEnergy_get(int Z)

◆ getMaterialParameters() [2/4]

void genfit::TGeoMaterialInterface::getMaterialParameters ( double &  density,
double &  Z,
double &  A,
double &  radiationLength,
double &  mEE 
)
virtual

Get material parameters in current material.

Implements genfit::AbsMaterialInterface.

◆ getMaterialParameters() [3/4]

void genfit::TGeoMaterialInterface::getMaterialParameters ( MaterialProperties parameters)
virtual

Implements genfit::AbsMaterialInterface.

Definition at line 74 of file TGeoMaterialInterface.cc.

74 {
75
76 TGeoMaterial* mat = gGeoManager->GetCurrentVolume()->GetMedium()->GetMaterial();
77
78 parameters.setMaterialProperties(mat->GetDensity(),
79 mat->GetZ(),
80 mat->GetA(),
81 mat->GetRadLen(),
83
84}

◆ getMaterialParameters() [4/4]

void genfit::TGeoMaterialInterface::getMaterialParameters ( MaterialProperties parameters)
virtual

◆ initTrack() [1/2]

bool genfit::TGeoMaterialInterface::initTrack ( double  posX,
double  posY,
double  posZ,
double  dirX,
double  dirY,
double  dirZ 
)
virtual

Initialize the navigator at given position and with given direction. Returns true if the volume changed.

Implements genfit::AbsMaterialInterface.

Definition at line 39 of file TGeoMaterialInterface.cc.

40 {
41 #ifdef DEBUG
42 std::cout << "TGeoMaterialInterface::initTrack. \n";
43 std::cout << "Pos "; TVector3(posX, posY, posZ).Print();
44 std::cout << "Dir "; TVector3(dirX, dirY, dirZ).Print();
45 #endif
46
47 // Move to the new point.
48 bool result = !gGeoManager->IsSameLocation(posX, posY, posZ, kTRUE);
49 // Set the intended direction.
50 gGeoManager->SetCurrentDirection(dirX, dirY, dirZ);
51 return result;
52}

◆ initTrack() [2/2]

bool genfit::TGeoMaterialInterface::initTrack ( double  posX,
double  posY,
double  posZ,
double  dirX,
double  dirY,
double  dirZ 
)
virtual

Initialize the navigator at given position and with given direction. Returns true if the volume changed.

Implements genfit::AbsMaterialInterface.


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