TimeOfFlight
1.32.00
|
TOFEstimators is the Marlin processor that calculates time-of-flight for charged particles that can be used for the mass reconstruction and particle identification.
The detailed description of the source code is available in the doxygen documentation.
If formulas are barely visible please switch to the light theme of the GitHub or check the main page in the doxygen documentation.
Time-of-flight is calculated only for PFOs with exactly one Track and exactly one Cluster. In other cases we write time-of-flight as 0.
TOFEstimators works with REC files but does not work with DST or mini-DST files as it requires individual hit information from the calorimeter.
You need to run InitDD4hep processor before running TOFEstimator. TOFEstimators is dependent on detector geometry to identify SET hits.
TOFEstimators has five steering parameters:
ReconstructedParticle
objects to analyze.ExtrapolateToEcal | default: true
If true
time-of-flight is calculated using ECal hits via one of the methods chosen in TofMethod steering parameter.
If false
time-of-flight is calculated using the last tracker hit. Which is most likely SET hit for the barrel region.
If the last tracker hit is the SET hit time-of-flight is set to the average time of two SET strips.
If the last tracker hit is not the SET hit time-of-flight is set to 0
(endcap case).
TofMethod | options: "closest", "frankAvg", "frankFit" | default: "closest"
If ExtrapolateToEcal is set to false
then this steering parameter is ignored.
If ExtrapolateToEcal is set to true
then it defines how to use ECal hits to calculate the time-of-flight to the ECal surface.
0.0
.0.0
.0.0
.for more illustrative explanations of the methods above see slides 10, 12, 13 from the following LCWS2021 talk).
TimeResolution | default: 0.0 (ps)
If ExtrapolateToEcal is set to true
it defines times resolution of individual ECal hits in ps.
If ExtrapolateToEcal is set to false
it defines times resolution of individual SET strips in ps.
MaxEcalLayer | default: 10
Defines number of inner ECal layers to use for the frankAvg and frankFit algorithms.
If frankAvg and frankFit are unused then this steering parameter is ignored.
This processor has single output parameter - time-of-flight.
./xml/steer.xml is a steering file example that runs three TOFEstimators processors: MyTofClosest0ps, MyTofSET10ps and MyTofFrankAvg50ps. They write output parameters in the PIDHandlers of the PandoraPFOs collection. Results are saved in the new output.slcio file using LCIOOutputProcessor.
To run this example one needs to setup iLCSoft environment. If the reader has a NAF account he/she can setup iLCSoft environment with:
source /cvmfs/ilc.desy.de/sw/x86_64_gcc82_centos7/v02-02-03/init_ilcsoft.sh
and run the example steering file with:
Marlin ./xml/steer.xml
Then one can look at the output.sclio file with, e.g.:
dumpevent output.slcio 1 | less
One can find output for PandoraPFOs which has new TOF algorithms attached...
collection name : PandoraPFOs parameters: --------------- print out of ReconstructedParticle collection --------------- parameter ParameterNames_MyTofClosest0ps [string]: timeOfFlight, parameter ParameterNames_MyTofFrankAvg50ps [string]: timeOfFlight, parameter ParameterNames_MyTofSET10ps [string]: timeOfFlight,
... and see final results for each individual PFO
------------ detailed PID info: --- algorithms : [id: 9] MyTofClosest0ps - params: timeOfFlight [id: 11] MyTofFrankAvg50ps - params: timeOfFlight [id: 10] MyTofSET10ps - params: timeOfFlight [particle] | PDG | likelihood | type | algoId | parameters : | | | | | [00000073] . . . | 0 | 0.0000e+00 | 000000 | 9 | [ timeOfFlight : 9.58e+00,] | 0 | 0.0000e+00 | 000000 | 10 | [ timeOfFlight : 0.00e+00,] | 0 | 0.0000e+00 | 000000 | 11 | [ timeOfFlight : 9.57e+00,]
After you run TOFEstimators and TrackLength processors you might want to run you analysis processor to e.g. calculate the mass of particles using time-of-flight information.
Here is the code example how to do that:
float YourAmazingAnalysisProcessor::getParameterFromPID(ReconstructedParticle* pfo, PIDHandler& pidHandler, std::string algorithmName, std::string parameterName){ int algorithmID = pidHandler.getAlgorithmID(algorithmName); const ParticleID& pfoPID = pidHandler.getParticleID(pfo, algorithmID); const std::vector<float>& parameters = pfoPID.getParameters(); int parIdx = pidHandler.getParameterIndex(algorithmID, parameterName); return parameters[parIdx]; } void YourAmazingAnalysisProcessor::processEvent(LCEvent* event){ LCCollection* pfos = event->getCollection("PandoraPFOs"); PIDHandler pidHandler(pfos); for(int i=0; i < pfos->getNumberOfElements(); ++i){ ReconstructedParticle* pfo = static_cast <ReconstructedParticle*> ( pfos->getElementAt(i) ); float momentum = getParameterFromPID(pfo, pidHandler, "MyTrackLengthProcessor", "momentumHMToEcal"); // in GeV float trackLength = getParameterFromPID(pfo, pidHandler, "MyTrackLengthProcessor", "trackLengthToEcal"); // in mm float tof = getParameterFromPID(pfo, pidHandler, "MyTofClosest0ps", "TimeOfFlight"); // in ns //calculate mass in GeV using relativistic momentum formula double mass = momentum * std::sqrt( std::pow(tof*CLHEP::c_light/trackLength, 2) - 1 ); } }