SND@LHC Software
Loading...
Searching...
No Matches
main.cc File Reference
#include <ConstField.h>
#include <Exception.h>
#include <FieldManager.h>
#include <KalmanFitterRefTrack.h>
#include <StateOnPlane.h>
#include <Track.h>
#include <TrackPoint.h>
#include <MaterialEffects.h>
#include <RKTrackRep.h>
#include <TGeoMaterialInterface.h>
#include <EventDisplay.h>
#include <HelixTrackModel.h>
#include <MeasurementCreator.h>
#include <TDatabasePDG.h>
#include <TEveManager.h>
#include <TGeoManager.h>
#include <TRandom.h>
#include <TVector3.h>
#include <vector>
#include <TMath.h>
Include dependency graph for main.cc:

Go to the source code of this file.

Functions

int main ()
 

Function Documentation

◆ main()

int main ( )

Definition at line 31 of file main.cc.

31 {
32
33 gRandom->SetSeed(14);
34
35 // init MeasurementCreator
36 genfit::MeasurementCreator measurementCreator;
37
38
39 // init geometry and mag. field
40 new TGeoManager("Geometry", "Geane geometry");
41 TGeoManager::Import("genfitGeom.root");
42 genfit::FieldManager::getInstance()->init(new genfit::ConstField(0.,0., 15.)); // 15 kGauss
44
45
46 // init event display
48
49
50 // init fitter
52
53 // main loop
54 for (unsigned int iEvent=0; iEvent<100; ++iEvent){
55
56 // true start values
57 TVector3 pos(0, 0, 0);
58 TVector3 mom(1.,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.));
62
63
64 // helix track model
65 const int pdg = 13; // particle pdg code
66 const double charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/(3.);
67 genfit::HelixTrackModel* helix = new genfit::HelixTrackModel(pos, mom, charge);
68 measurementCreator.setTrackModel(helix);
69
70
71 unsigned int nMeasurements = gRandom->Uniform(5, 15);
72
73
74 // smeared start values
75 const bool smearPosMom = true; // init the Reps with smeared pos and mom
76 const double posSmear = 0.1; // cm
77 const double momSmear = 3. /180.*TMath::Pi(); // rad
78 const double momMagSmear = 0.1; // relative
79
80 TVector3 posM(pos);
81 TVector3 momM(mom);
82 if (smearPosMom) {
83 posM.SetX(gRandom->Gaus(posM.X(),posSmear));
84 posM.SetY(gRandom->Gaus(posM.Y(),posSmear));
85 posM.SetZ(gRandom->Gaus(posM.Z(),posSmear));
86
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()));
90 }
91 // approximate covariance
92 TMatrixDSym covM(6);
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);
98
99
100 // trackrep
102
103 // smeared start state
104 genfit::MeasuredStateOnPlane stateSmeared(rep);
105 rep->setPosMomCov(stateSmeared, posM, momM, covM);
106
107
108 // create track
109 TVectorD seedState(6);
110 TMatrixDSym seedCov(6);
111 rep->get6DStateCov(stateSmeared, seedState, seedCov);
112 genfit::Track fitTrack(rep, seedState, seedCov);
113
114
115 // create random measurement types
116 std::vector<genfit::eMeasurementType> measurementTypes;
117 for (unsigned int i = 0; i < nMeasurements; ++i)
118 measurementTypes.push_back(genfit::eMeasurementType(gRandom->Uniform(8)));
119
120
121 // create smeared measurements and add to track
122 try{
123 for (unsigned int i=0; i<measurementTypes.size(); ++i){
124 std::vector<genfit::AbsMeasurement*> measurements = measurementCreator.create(measurementTypes[i], i*5.);
125 fitTrack.insertPoint(new genfit::TrackPoint(measurements, &fitTrack));
126 }
127 }
128 catch(genfit::Exception& e){
129 std::cerr<<"Exception, next track"<<std::endl;
130 std::cerr << e.what();
131 continue;
132 }
133
134 //check
135 assert(fitTrack.checkConsistency());
136
137 // do the fit
138 fitter->processTrack(&fitTrack);
139
140 //check
141 assert(fitTrack.checkConsistency());
142
143
144 if (iEvent < 1000) {
145 // add track to event display
146 display->addEvent(&fitTrack);
147 }
148
149
150
151 }// end loop over events
152
153 delete fitter;
154
155 // open event display
156 display->open();
157
158}
Abstract base class for Kalman fitter and derived fitting algorithms.
Abstract base class for a track representation.
Definition AbsTrackRep.h:66
virtual void setPosMomCov(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const =0
Set position, momentum and covariance of state.
virtual void get6DStateCov(const MeasuredStateOnPlane &state, TVectorD &stateVec, TMatrixDSym &cov) const
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
Constant Magnetic field.
Definition ConstField.h:37
Event display designed to run with Genfit.
void addEvent(std::vector< genfit::Track * > &tracks)
Add new event.
static EventDisplay * getInstance()
void open()
Open the event display.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition Exception.h:48
virtual const char * what() const
Standard error message handling for exceptions. use like "std::cerr << e.what();".
Definition Exception.cc:51
void init(AbsBField *b)
set the magnetic field here. Magnetic field classes must be derived from AbsBField.
static FieldManager * getInstance()
Get singleton instance.
Helix track model for testing purposes.
Kalman filter implementation with linearization around a reference track.
static MaterialEffects * getInstance()
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.
std::vector< genfit::AbsMeasurement * > create(eMeasurementType, double tracklength, bool &outlier, int &lr)
void setTrackModel(const HelixTrackModel *model)
Takes ownership!
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u', v', u, v)
Definition RKTrackRep.h:70
AbsMaterialInterface implementation for use with ROOT's TGeoManager.
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition TrackPoint.h:50
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Definition Track.h:71
int i
Definition ShipAna.py:86