1 #include <ConstField.h>
3 #include <FieldManager.h>
4 #include <KalmanFitterRefTrack.h>
5 #include <StateOnPlane.h>
7 #include <TrackPoint.h>
9 #include <MaterialEffects.h>
10 #include <RKTrackRep.h>
11 #include <TGeoMaterialInterface.h>
13 #include <EventDisplay.h>
15 #include <HelixTrackModel.h>
16 #include <MeasurementCreator.h>
18 #include <TDatabasePDG.h>
19 #include <TEveManager.h>
20 #include <TGeoManager.h>
25 #include "TDatabasePDG.h"
40 new TGeoManager(
"Geometry",
"Geane geometry");
41 TGeoManager::Import(
"genfitGeom.root");
54 for (
unsigned int iEvent=0; iEvent<100; ++iEvent){
57 TVector3 pos(0, 0, 0);
59 mom.SetPhi(gRandom->Uniform(0.,2*TMath::Pi()));
60 mom.SetTheta(gRandom->Uniform(0.4*TMath::Pi(),0.6*TMath::Pi()));
61 mom.SetMag(gRandom->Uniform(0.2, 1.));
66 const double charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/(3.);
71 unsigned int nMeasurements = gRandom->Uniform(5, 15);
75 const bool smearPosMom =
true;
76 const double posSmear = 0.1;
77 const double momSmear = 3. /180.*TMath::Pi();
78 const double momMagSmear = 0.1;
83 posM.SetX(gRandom->Gaus(posM.X(),posSmear));
84 posM.SetY(gRandom->Gaus(posM.Y(),posSmear));
85 posM.SetZ(gRandom->Gaus(posM.Z(),posSmear));
87 momM.SetPhi(gRandom->Gaus(mom.Phi(),momSmear));
88 momM.SetTheta(gRandom->Gaus(mom.Theta(),momSmear));
89 momM.SetMag(gRandom->Gaus(mom.Mag(), momMagSmear*mom.Mag()));
93 double resolution = 0.01;
94 for (
int i = 0; i < 3; ++i)
95 covM(i,i) = resolution*resolution;
96 for (
int i = 3; i < 6; ++i)
97 covM(i,i) = pow(resolution / nMeasurements /
sqrt(3), 2);
105 stateSmeared.setPosMomCov(posM, momM, covM);
109 TVectorD seedState(6);
110 TMatrixDSym seedCov(6);
111 stateSmeared.get6DStateCov(seedState, seedCov);
116 std::vector<genfit::eMeasurementType> measurementTypes;
117 for (
unsigned int i = 0; i < nMeasurements; ++i)
118 measurementTypes.push_back(genfit::eMeasurementType(gRandom->Uniform(8)));
123 for (
unsigned int i=0; i<measurementTypes.size(); ++i){
124 std::vector<genfit::AbsMeasurement*> measurements = measurementCreator.create(measurementTypes[i], i*5.);
129 std::cerr<<
"Exception, next track"<<std::endl;
130 std::cerr << e.what();
135 fitTrack.checkConsistency();
138 fitter->processTrack(&fitTrack);
141 fitTrack.checkConsistency();
146 display->addEvent(&fitTrack);
Abstract base class for Kalman fitter and derived fitting algorithms.
Abstract base class for a track representation.
Event display designed to run with Genfit.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
static FieldManager * getInstance()
Get singleton instance.
void init(AbsBField *b)
set the magnetic field here. Magnetic field classes must be derived from AbsBField.
Helix track model for testing purposes.
Kalman filter implementation with linearization around a reference track.
void init(AbsMaterialInterface *matIfc)
set the material interface here. Material interface classes must be derived from AbsMaterialInterface...
#StateOnPlane with additional covariance matrix.
Create different measurement types along a HelixTrackModel for testing purposes.
void setTrackModel(const HelixTrackModel *model)
Takes ownership!
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u', v', u, v)
AbsMaterialInterface implementation for use with ROOT's TGeoManager.
Object containing AbsMeasurement and AbsFitterInfo objects.
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
double sqrt(double a)
sqrt for double
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
int main(int argc, char **argv)
Run all tests.