]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Possibility to define the magnetic field in the reconstruction (Yu.Belikov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 17 Apr 2001 08:07:14 +0000 (08:07 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 17 Apr 2001 08:07:14 +0000 (08:07 +0000)
12 files changed:
ITS/AliITSrecoV2.h
ITS/AliITStestV2.C
ITS/AliITStrackV2.cxx
ITS/AliITStrackV2.h
STEER/AliKalmanTrack.cxx
STEER/AliKalmanTrack.h
TPC/AliTPCComparison.C
TPC/AliTPCclusterer.cxx
TPC/AliTPCtest.C
TPC/AliTPCtrack.cxx
TPC/AliTPCtrack.h
TPC/AliTPCtracker.cxx

index a171b9c3483098eb92f75c25144ab4f1c4364165..93603a7389ae8dafebabcff34300dcd0bb0819a7 100644 (file)
@@ -9,7 +9,6 @@
 //       Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
 //-------------------------------------------------------------------------
 #include <Rtypes.h>
-#include <iostream.h>
 
 //namespace AliITSreco {    
    const Int_t kMaxClusterPerLayer=3500*10;
@@ -31,8 +30,6 @@
 
    const Double_t kSigmaYV=0.005e+0;
    const Double_t kSigmaZV=0.010e+0;
-
-   const Double_t kConvConst=100/0.299792458/0.2; 
 //}
 
 //using namespace AliITSreco;
index 197a120545b3b988fa20ee881a53ce6e4be24475..91160ce7dc7e76224bd7369a4174dc7f08948975 100644 (file)
@@ -1,6 +1,20 @@
 Int_t AliITStestV2() {
    Int_t rc=0;
 
+   if (gAlice) {delete gAlice; gAlice=0;}
+   TFile *in=TFile::Open("galice.root");
+   if (!in->IsOpen()) {
+      cerr<<"Can't open galice.root !\n"; 
+      return 1;
+   }
+   if (!(gAlice=(AliRun*)in->Get("gAlice"))) {
+      cerr<<"Can't find gAlice !\n";
+      return 2;
+   }
+   AliKalmanTrack::SetConvConst(100/0.299792458/0.2/gAlice->Field()->Factor());
+   delete gAlice; gAlice=0;
+   in->Close();
+
    gROOT->LoadMacro("$(ALICE_ROOT)/ITS/AliITSFindClustersV2.C");
    if (rc=AliITSFindClustersV2()) return rc;
 
index 39024ac1f26a4491cfe338ef8cb82ffe462afcad..193ab605da02a4ae6e6cc5cdf28916251d742bec 100644 (file)
@@ -35,7 +35,7 @@ const Int_t kWARN=1;
 //____________________________________________________________________________
 AliITStrackV2::AliITStrackV2(const AliTPCtrack& t) throw (const Char_t *) {
   //------------------------------------------------------------------
-  //Convertion TPC track -> ITS track
+  //Conversion TPC track -> ITS track
   //------------------------------------------------------------------
   SetLabel(t.GetLabel());
   SetChi2(0.);
@@ -45,16 +45,16 @@ AliITStrackV2::AliITStrackV2(const AliTPCtrack& t) throw (const Char_t *) {
   if      (fAlpha < -TMath::Pi()) fAlpha += 2*TMath::Pi();
   else if (fAlpha >= TMath::Pi()) fAlpha -= 2*TMath::Pi();
 
-  //Convertion of the track parameters
+  //Conversion of the track parameters
   Double_t x,p[5]; t.GetExternalParameters(x,p);
-  fX=x;    x=kConversionConstant;
+  fX=x;    x=GetConvConst();
   fP0=p[0]; 
   fP1=p[1]; 
   fP2=p[2];
   fP3=p[3];
   fP4=p[4]/x; 
 
-  //Convertion of the covariance matrix
+  //Conversion of the covariance matrix
   Double_t c[15]; t.GetExternalCovariance(c);
 
   fC00=c[0 ];
@@ -106,7 +106,7 @@ void AliITStrackV2::GetExternalCovariance(Double_t cc[15]) const {
   // This function returns an external representation of the covriance matrix.
   //   (See comments in AliTPCtrack.h about external track representation)
   //-------------------------------------------------------------------------
-  Double_t a=kConvConst;
+  Double_t a=GetConvConst();
 
   cc[0 ]=fC00;
   cc[1 ]=fC10;   cc[2 ]=fC11;
index 2defb7acdd70606a71978f93eb0dec9afcde8fc4..136507365047c33c41fe0098b130ccc6eb295bb4 100644 (file)
@@ -53,7 +53,7 @@ public:
   Double_t GetZ()    const {return fP1;}
   Double_t GetSnp()  const {return fP2;}
   Double_t GetTgl()  const {return fP3;}
-  Double_t Get1Pt()  const {return fP4*kConvConst;}
+  Double_t Get1Pt()  const {return fP4*GetConvConst();}
 
 
   Double_t GetD() const;
index e5c43794c8820fcd8f78247f8dbbd407f645f149..c82f383ea5cd6dbf33d1fe00bcfd1968fd55222f 100644 (file)
@@ -23,3 +23,5 @@
 
 ClassImp(AliKalmanTrack)
 
+Double_t AliKalmanTrack::fConvConst;
+
index a4da4045b59f3e941365f6ed0304547499d85d6e..b23c84534bdd92bbe217213a142b25dc031d12dc 100644 (file)
@@ -16,7 +16,7 @@ class AliCluster;
 
 class AliKalmanTrack : public TObject {
 public:
-  AliKalmanTrack() {fLab=-3141593; fChi2=0; fN=0;}
+  AliKalmanTrack() { fLab=-3141593; fChi2=0; fN=0; }
   AliKalmanTrack(const AliKalmanTrack &t) {fLab=t.fLab;fChi2=t.fChi2;fN=t.fN;}
   virtual ~AliKalmanTrack(){};
   void SetLabel(Int_t lab) {fLab=lab;}
@@ -40,6 +40,9 @@ public:
   Int_t PropagateTo(Double_t xr,Double_t x0,Double_t rho,Double_t pm)=0;
   virtual Int_t Update(const AliCluster* c, Double_t chi2, UInt_t i)=0;
 
+  static void SetConvConst(Double_t cc) {fConvConst=cc;}
+  Double_t GetConvConst() const {return fConvConst;}
+
 protected:
   void SetChi2(Double_t chi2) {fChi2=chi2;} 
   void SetNumberOfClusters(Int_t n) {fN=n;} 
@@ -49,6 +52,8 @@ private:
   Double_t fChi2;         // total chi2 value for this track
   Int_t fN;               // number of associated clusters
 
+  static Double_t fConvConst; //conversion constant cm -> GeV/c
+
   ClassDef(AliKalmanTrack,1)    // Reconstructed track
 };
 
index c595a5582d697bf3b45a67d2ae9978ea1636335d..0e8f92f636590d1c423bd731934ac18f800c9217 100644 (file)
@@ -11,7 +11,8 @@ struct GoodTrack {
 Int_t good_tracks(GoodTrack *gt, Int_t max);
 
 Int_t AliTPCComparison() {
-  Int_t i;
+   cerr<<"Doing comparison...\n";
+   Int_t i;
    gBenchmark->Start("AliTPCComparison");
 
    TFile *cf=TFile::Open("AliTPCclusters.root");
@@ -92,7 +93,6 @@ Int_t AliTPCComparison() {
    }
    cerr<<"Number of good tracks : "<<ngood<<endl;
 
-   cerr<<"Doing comparison...\n";
    TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
    TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
    TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
@@ -109,7 +109,7 @@ Int_t AliTPCComparison() {
    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
 
    TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
-   TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,200.);
+   TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,250.);
    hep->SetMarkerStyle(8);
    hep->SetMarkerSize(0.4);
 
@@ -257,9 +257,9 @@ Int_t AliTPCComparison() {
    hg->SetXTitle("Pt (GeV/c)");
    hg->Draw();
 
-   TLine *line1 = new TLine(0,1.0,7,1.0); line1->SetLineStyle(4);
+   TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
    line1->Draw("same");
-   TLine *line2 = new TLine(0,0.9,7,0.9); line2->SetLineStyle(4);
+   TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
    line2->Draw("same");
 
    hf->SetFillColor(1);
@@ -300,7 +300,8 @@ Int_t AliTPCComparison() {
 Int_t good_tracks(GoodTrack *gt, Int_t max) {
    Int_t nt=0;
 
-   TFile *file=TFile::Open("rfio:galice.root");
+   //TFile *file=TFile::Open("rfio:galice.root");
+   TFile *file=TFile::Open("galice.root");
    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
 
    if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
@@ -379,7 +380,7 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
           digp->AdjustSectorRow(digits->GetID(),sec,row);
           cerr<<sec<<' '<<row<<"                                     \r";
           digits->First();
-          while (digits->Next()) {
+          do { //Many thanks to J.Chudoba who noticed this
               Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
               Short_t dig = digits->GetDigit(it,ip);
               Int_t idx0=digits->GetTrackID(it,ip,0); 
@@ -388,7 +389,7 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
               if (idx0>=0 && dig>=zero) count[idx0]+=1;
               if (idx1>=0 && dig>=zero) count[idx1]+=1;
               if (idx2>=0 && dig>=zero) count[idx2]+=1;
-          }
+          } while (digits->Next());
           for (Int_t j=0; j<np; j++) {
               if (count[j]>1) {
                  if (sec>=digp->GetNInnerSector())
index b01dc4d8bf6ec56bee9b4a893f2bdcfd185ab555..3980b9f22fb3b5f16a0a62a90267808ede03abe2 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.3  2000/10/05 16:14:01  kowal2
+Forward declarations.
+
 Revision 1.2  2000/06/30 12:07:50  kowal2
 Updated from the TPC-PreRelease branch
 
@@ -242,8 +245,8 @@ void AliTPCclusterer::Digits2Clusters(const AliTPCParam *par, TFile *of)
     carray.StoreRow(sec,row);
     carray.ClearRow(sec,row);
 
-    cerr<<"sector, row, compressed digits, clusters: "
-    <<sec<<' '<<row<<' '<<digarr.GetSize()<<' '<<ncl<<"                  \r";
+    //cerr<<"sector, row, compressed digits, clusters: "
+    //<<sec<<' '<<row<<' '<<digarr.GetSize()<<' '<<ncl<<"                  \r";
 
     nclusters+=ncl;
 
index 975625e688aaa610a04ceb6ef6fd014b08e12cd4..f21306c82246fbf6fdf995991e71e0404d7b6304 100644 (file)
@@ -5,6 +5,8 @@ Int_t AliTPCtest() {
    gROOT->LoadMacro("$(ALICE_ROOT)/macros/grun.C");
    grun();
 
+   AliKalmanTrack::SetConvConst(100/0.299792458/0.2/gAlice->Field()->Factor());
+
    Int_t ver=gAlice->GetDetector("TPC")->IsVersion();
    delete gAlice; gAlice=0;
 
index dedeba16cc3f47a18207bbd16f4fd430006e48fa..c3a67c059b2da18e603e79ad32538a6dae68cbb2 100644 (file)
@@ -92,7 +92,7 @@ void AliTPCtrack::GetExternalCovariance(Double_t cc[15]) const {
   // This function returns an external representation of the covriance matrix.
   //   (See comments in AliTPCtrack.h about external track representation)
   //-------------------------------------------------------------------------
-  Double_t a=kConversionConstant;
+  Double_t a=GetConvConst();
 
   Double_t c22=fX*fX*fC33-2*fX*fC32+fC22;
   Double_t c42=fX*fC43-fC42;
index 36a2ea6c250af61deccf01b945f6806a79255614..30b30b75662e2ae316c361a0f5139d62ef7dd5f2 100644 (file)
@@ -26,8 +26,6 @@
 #include <AliKalmanTrack.h>
 #include <TMath.h>
 
-const Double_t kConversionConstant=100/0.299792458/0.2; 
-
 class AliTPCClustersArray;
 class AliTPCcluster;
 
@@ -51,7 +49,7 @@ public:
   Double_t GetY()   const {return fP0;}
   Double_t GetZ()   const {return fP1;}
   Double_t GetSnp() const {return fX*fP3 - fP2;}             
-  Double_t Get1Pt() const {return fP3*kConversionConstant;}             
+  Double_t Get1Pt() const {return fP3*GetConvConst();}             
   Double_t GetTgl() const {return fP4;}
 
   Double_t GetSigmaY2() const {return fC00;}
index 5cc4601b4e3f64a0731374d6c8bac8a938cd6211..e7976c58833551421f07689631e470abe60dc761 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.6  2001/03/13 14:25:47  hristov
+New design of tracking classes (Yu.Belikov)
+
 Revision 1.5  2000/12/20 07:51:59  kowal2
 Changes suggested by Alessandra and Paolo to avoid overlapped
 data fields in encapsulated classes.
@@ -431,7 +434,7 @@ Int_t AliTPCtracker::Clusters2Tracks(const AliTPCParam *par, TFile *of) {
              iotrack=pt;
              tracktree.Fill();
              t.UseClusters(&carray,nc);
-             cerr<<found++<<'\r';
+             found++;
           }
        }
     }