kalmanalignment is hosted by Hepforge, IPPP Durham

/home/mweber/Dokumente/CMS/Alignment/kalmanalignment/trunk/simulation/simulation.cxx File Reference

This is the main file for simulation and track reconstruction. More...

#include <iostream>
#include <iomanip>
#include <string>
#include <vector>
#include "TROOT.h"
#include "TRint.h"
#include "TRandom.h"
#include "TBenchmark.h"
#include "TCanvas.h"
#include "TH1.h"
#include "TF1.h"
#include "TEnv.h"
#include "TMath.h"
#include "TStyle.h"
#include "TVectorD.h"
#include "TMatrixDSym.h"
#include "TFile.h"
#include "TTree.h"
#include "TArray.h"
#include "CLHEP/Matrix/Vector.h"
#include "CLHEP/Matrix/Matrix.h"
#include "CLHEP/Matrix/SymMatrix.h"
#include "Utilities.h"
#include "Constants.h"
#include "Det.h"
#include "Detector.h"
#include "Track.h"
#include "ThreeDModel.h"
#include "Hit.h"
#include "KalmanFilterTracking.h"
#include "Visualization.h"
#include "AlignEvent.h"
#include "KalmanFilterAlignmentInputProvider.h"

Enumerations

enum  DetectorResolutionFunction { kDelta, kUniform, kGaus }
 

How to simulate detector resolution.

More...

Functions

Track generate_cosmic_track ()
 simulate gaussian resolution
Track generate_uniform_track (double momentum=4)
 Generate a track, uniformely distributed over the whole detector area.
void check_measurement_model (const Track &track, const vector< Hit * > hitvector)
 Check the measurement model.
double ms_angle (double beta, int z, double p, double x, double xi)
 Calculate the projected mean angular deflection due to multiple scattering of a particle in material.
void simulate_multiple_scattering (Track &t, const Det &adet, const Hit *hit)
 Simulation of multiple scattering.
void simulate_hits (Track &t, DetectorResolutionFunction smear=kUniform, bool multiple_scattering=false)
 Simulation step: Simulate hits in real (misaligned) detector.
void check_track_parameter_transformation (Track t)
 Test transformation between global and local track parameters.
void check_track_covariance_transformation ()
 Test transformation between global and local track covariance.
void check_track_linearization ()
void fill_tracking_pulls (const Track &reco, const Track &truth)
 Compute tracking pulls, difference to truth dividided by tracking errors.
void fill_track_chi2 (const Track &reco)
 Make a histogram of track chi2, chi2/ndof and chi2-Probability.
void subtest ()
int main (int argc, char *argv[])

Variables

Visualization gVis
 Global Visualization object.
double gMinimumEnergy
 Global configuration for minimum energy (e.g. lead shield simulation)

Detailed Description

This is the main file for simulation and track reconstruction.


Enumeration Type Documentation

How to simulate detector resolution.

Enumerator:
kUniform 

perfect detector resolution

kGaus 

simulate uniform resolution


Function Documentation

void check_measurement_model ( const Track track,
const vector< Hit * >  hitvector 
)

Check the measurement model.

This function compares the predicted hit (intersection of track with the nominal detector position) with the calculation from the measurement model. If the true track parameters and the true misalignment is known, the measurement model allows to compute the true position from the predicted position.

The values reported should be identical. However, this is only true if there has been no simulation of detector resolution and no multiple scattering effects are being computed. If those are incorporated, differences are expected.

void check_track_linearization ( )

< plot nPoints points

< plot in [-range,range) around current track parameters

create canvas for plotting

This is where the track linearization takes place

loop over all detectors

get detector

compute track parameters at detector location

loop over all four parameters and fill extrapolation histograms

create histograms

< add small step to avoid strange binning effect (center of bin)

fill histograms

void fill_tracking_pulls ( const Track reco,
const Track truth 
)

Compute tracking pulls, difference to truth dividided by tracking errors.

This function uses the known truth, i.e. the true track parameters, and compares the reconstructed ones to them, dividing by the track errors. This "pull" distribution should be, in case of correct parameters and error estimates, have mean zero and RMS one. With gaussian input errors (i.e. gaussian hit errors), it should follow a gaussian distribution. However, hits are distributed according to uniform distribution, therefore the shape is not gaussian (but mean should still be zero and RMS one).

Track generate_cosmic_track ( )

simulate gaussian resolution

Generate a track, distributed according to the cosmic muon spectrum over the whole detector area.

The muon momentum distribution is taken from http://arxiv.org/pdf/hep-ph/0103322.pdf and has been fitted by a third order polynomial. Generation is only done in a cone with an opening angle of 20 degrees, according to the Caprice measurements.

Track generate_uniform_track ( double  momentum = 4)

Generate a track, uniformely distributed over the whole detector area.

Parameters:
momentumis the track momentum. It default value is 4 GeV, which is the mean energy of the cosmic muon momentum distribution.
int main ( int  argc,
char *  argv[] 
)

reconstruct track using nominal coordinate system

double ms_angle ( double  beta,
int  z,
double  p,
double  x,
double  xi 
)

Calculate the projected mean angular deflection due to multiple scattering of a particle in material.

Projected means that the angle is not the angle in space, but rather a one-dimensional projection on one axis of a plane that is perpendicular to the particle direction. To obtaint the spatial angle, one has to use twice the returned value.

Parameters:
betaVelocity/c
zCharge of particle)
pMomentum, must be given in MeV
xMatter thickness in cm
xiRadiation length of material in cm
void simulate_hits ( Track t,
DetectorResolutionFunction  smear = kUniform,
bool  multiple_scattering = false 
)

Simulation step: Simulate hits in real (misaligned) detector.

Compute hits in all detectors for the given track. Tracks are smeared with the detector resolution if asked for.

void simulate_multiple_scattering ( Track t,
const Det adet,
const Hit hit 
)

Simulation of multiple scattering.

The track is scattered at the given detector. The offset through scattering in the material is neglected.

Parameters:
tThe track before (and after) multiple scattering. The track parameters are modified according to the new direction.
adetDetector in which the multiple scattering occurs.
hitMust be the hit before applying any other smearing, because it is taken as the starting point for the resulting track!

The routine estimates the radiation length according to the path length of the particle track through the silicon detector.

void subtest ( )

Test HepMatrix.sub()