Index: trunk/MagicSoft/Mars/Changelog
===================================================================
--- trunk/MagicSoft/Mars/Changelog	(revision 1236)
+++ trunk/MagicSoft/Mars/Changelog	(revision 1237)
@@ -1,3 +1,11 @@
                                                                   -*-*- END -*-*-
+
+ 2002/03/07: Thomas Bretz
+
+   * manalysis/MHillasSrc.[h,cc]:
+     - added fHeadTail
+     - changed version number to 2.
+
+
 
  2002/03/07: Thomas Bretz
Index: trunk/MagicSoft/Mars/manalysis/MHillasSrc.cc
===================================================================
--- trunk/MagicSoft/Mars/manalysis/MHillasSrc.cc	(revision 1236)
+++ trunk/MagicSoft/Mars/manalysis/MHillasSrc.cc	(revision 1237)
@@ -69,11 +69,14 @@
     const Double_t sy = my - fSrcPos->GetY();      // [mm]
 
-    const Double_t tand = tan(hillas->GetDelta()); // [1]
+    const Double_t sd = sin(hillas->GetDelta());   // [1]
+    const Double_t cd = cos(hillas->GetDelta());   // [1]
 
-    fDist   = sqrt(sx*sx + sy*sy);                 // [mm]
-    fAlpha  = atan2(sy - tand*sx, sx + tand*sy);   // [rad]
-    fAlpha *= kRad2Deg;                            // [deg]
 
-    fHillas = hillas;
+    fHeadTail = cd*sx + sd*sy;                     // [mm]
+    fDist     = sqrt(sx*sx + sy*sy);               // [mm]
+    fAlpha    = atan2(cd*sy - sd*sx, fHeadTail);   // [rad]
+    fAlpha   *= kRad2Deg;                          // [deg]
+
+    fHillas   = hillas;
 
     SetReadyToSave();
@@ -84,6 +87,7 @@
     *fLog << all;
     *fLog << "Source dependant Image Parameters (" << GetName() << ")" << endl;
-    *fLog << " - Dist  = " << fDist  << " mm"  << endl;
-    *fLog << " - Alpha = " << fAlpha << " deg" << endl;
+    *fLog << " - Dist     = " << fDist     << " mm"  << endl;
+    *fLog << " - Alpha    = " << fAlpha    << " deg" << endl;
+    *fLog << " - HeadTail = " << fHeadTail << " mm"  << endl;
 }
 
@@ -96,4 +100,5 @@
     fin >> fAlpha;
     fin >> fDist;
+    fin >> fHeadTail;
 }
 
Index: trunk/MagicSoft/Mars/manalysis/MHillasSrc.h
===================================================================
--- trunk/MagicSoft/Mars/manalysis/MHillasSrc.h	(revision 1236)
+++ trunk/MagicSoft/Mars/manalysis/MHillasSrc.h	(revision 1237)
@@ -14,6 +14,7 @@
     const MSrcPosCam *fSrcPos; //! Source position in the camere
 
-    Float_t fAlpha;  // [deg]  angle of major axis with vector to src
-    Float_t fDist;   // [mm]   distance from src to center of ellipse
+    Float_t fAlpha;    // [deg]  angle of major axis with vector to src
+    Float_t fDist;     // [mm]   distance from src to center of ellipse
+    Float_t fHeadTail; // [mm]
 
 public:
@@ -23,12 +24,13 @@
     const MSrcPosCam *GetSrcPos() const   { return fSrcPos; }
 
-    Float_t GetLength() const  { return fHillas->GetLength(); }
-    Float_t GetWidth()  const  { return fHillas->GetWidth(); }
-    Float_t GetDelta()  const  { return fHillas->GetDelta(); }
-    Float_t GetSize()   const  { return fHillas->GetSize(); }
-    Float_t GetMeanX()  const  { return fHillas->GetMeanX(); }
-    Float_t GetMeanY()  const  { return fHillas->GetMeanY(); }
-    Float_t GetAlpha()  const  { return fAlpha; }
-    Float_t GetDist()   const  { return fDist; }
+    Float_t GetLength()   const  { return fHillas->GetLength(); }
+    Float_t GetWidth()    const  { return fHillas->GetWidth(); }
+    Float_t GetDelta()    const  { return fHillas->GetDelta(); }
+    Float_t GetSize()     const  { return fHillas->GetSize(); }
+    Float_t GetMeanX()    const  { return fHillas->GetMeanX(); }
+    Float_t GetMeanY()    const  { return fHillas->GetMeanY(); }
+    Float_t GetAlpha()    const  { return fAlpha; }
+    Float_t GetDist()     const  { return fDist; }
+    Float_t GetHeadTail() const  { return fHeadTail; }
 
     void Print(Option_t *opt=NULL) const;
@@ -39,6 +41,11 @@
     //virtual void AsciiWrite(ofstream &fout) const;
 
-    ClassDef(MHillasSrc, 1) // Container to hold source position dependant parameters
+    ClassDef(MHillasSrc, 2) // Container to hold source position dependant parameters
 };
 
+/*
+ Version 1: fAlpha, fDist
+ Version 2: fAlpha, fDist, fSign
+ */
+
 #endif
