]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliAstrolab.h
For Pythia with tune don't switch off MI in ConfigHeavyFlavor
[u/mrichter/AliRoot.git] / RALICE / AliAstrolab.h
1 #ifndef ALIASTROLAB_H
2 #define ALIASTROLAB_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 // $Id$
8
9 #include <math.h>
10
11 #include "TTask.h"
12 #include "TString.h"
13 #include "TRotMatrix.h"
14 #include "TObjArray.h"
15 #include "TArrayI.h"
16 #include "TCanvas.h"
17 #include "TEllipse.h"
18 #include "TLine.h"
19 #include "TText.h"
20 #include "TH2.h"
21 #include "TMarker.h"
22
23 #include "AliTimestamp.h"
24 #include "AliPosition.h"
25 #include "AliSignal.h"
26  
27 class AliAstrolab : public TTask,public AliTimestamp
28 {
29  public:
30   AliAstrolab(const char* name="AliAstrolab",const char* title="Generic lab"); // Constructor
31   virtual ~AliAstrolab();                                      // Destructor
32   AliAstrolab(const AliAstrolab& t);                           // Copy constructor
33   virtual TObject* Clone(const char* name="") const;           // Make a deep copy and provide its pointer
34   void Data(Int_t mode=1,TString u="deg");                     // Lab info in angular units u
35   void SetLabPosition(Ali3Vector& r);                          // Set lab position in terrestrial frame
36   void SetLabPosition(Double_t l,Double_t b,TString u="deg");  // Set lab terrestrial position
37   AliPosition GetLabPosition() const;                          // Provide the lab terrestrial position 
38   void GetLabPosition(Double_t& l,Double_t& b,TString u="deg") const;// Provide the lab terrestrial position
39   using AliTimestamp::GetLT;
40   Double_t GetLT();  // Provide Local Time (LT) in fractional hours
41   using AliTimestamp::GetLMST;
42   Double_t GetLMST(); // Provide Local Mean Sidereal Time (LMST) in fractional hours
43   using AliTimestamp::GetLAST;
44   Double_t GetLAST(); // Provide Local Apparent Sidereal Time (LAST) in fractional hours
45   using AliTimestamp::SetLT;
46   void SetLT(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns=0,Int_t ps=0); // Set specified LT
47   void SetLT(Int_t y,Int_t d,Int_t s,Int_t ns=0,Int_t ps=0); // Set LT based on elapsed days, secs etc...
48   Double_t ConvertAngle(Double_t a,TString in,TString out) const;       // Angular format conversions
49   void PrintAngle(Double_t a,TString in,TString out,Int_t ndig=1) const;// Print angle in various formats
50   void SetSignal(Ali3Vector* r,TString frame,TString mode,AliTimestamp* ts,Int_t jref=0,TString name=""); // Store a generic signal
51   void SetSignal(Double_t a,Double_t d,TString s,Double_t e,TString mode,Int_t jref=0,TString name="");   // Store RA, decl. and time
52   void SetSignal(Double_t a,Double_t d,TString mode,AliTimestamp* ts,Int_t jref=0,TString name="");       // Store RA, decl. and time
53   AliSignal* GetSignal(Ali3Vector& r,TString frame,TString mode,AliTimestamp* ts,Int_t jref=0);// Provide stored signal data
54   AliSignal* GetSignal(Ali3Vector& r,TString frame,TString mode,AliTimestamp* ts,TString name);// Provide stored signal data
55   AliSignal* GetSignal(Double_t& a,Double_t& d,TString mode,AliTimestamp* ts,Int_t jref=0);    // Provide corrected RA and decl.
56   AliSignal* GetSignal(Double_t& a,Double_t& d,TString mode,AliTimestamp* ts,TString name);    // Provide corrected RA and decl.
57   AliSignal* GetSignal(Double_t& a,Double_t& d,TString s,Double_t e,TString mode,Int_t jref=0);// Provide corrected RA and decl.
58   AliSignal* GetSignal(Double_t& a,Double_t& d,TString s,Double_t e,TString mode,TString name);// Provide corrected RA and decl.
59   AliSignal* GetSignal(Int_t jref=0);                // Provide pointer to a stored signal object
60   AliSignal* GetSignal(TString name);                // Provide pointer to a stored signal object
61   void RemoveRefSignal(Int_t j,Int_t compress);      // Remove a stored reference signal object
62   void RemoveRefSignal(TString name,Int_t compress); // Remove a stored reference signal object
63   void PrintSignal(TString frame,TString mode,AliTimestamp* ts,Int_t ndig,Int_t jref=0); // Print stored signal data
64   void PrintSignal(TString frame,TString mode,AliTimestamp* ts,Int_t ndig,TString name); // Print stored signal data
65   void ListSignals(TString frame,TString mode,Int_t ndig=1); // List all stored signals
66   Int_t GetSignalIndex(TString name); // Provide storage index of the signal with the specified name
67   Double_t GetHourAngle(TString mode,AliTimestamp* ts,Int_t jref=0);// Provide the Local Hour Angle in degrees
68   void SetLocalFrame(Double_t t1,Double_t p1,Double_t t2,Double_t p2,Double_t t3,Double_t p3); // Define local coordinate frame
69   using AliTimestamp::GetDifference;
70   Double_t GetDifference(Int_t jref,TString au,Double_t& dt,TString tu,Int_t mode=1,Int_t* ia=0,Int_t* it=0); // Provide space and time difference
71   Double_t GetDifference(TString name,TString au,Double_t& dt,TString tu,Int_t mode=1);// Provide space and time difference
72   TArrayI* MatchRefSignal(Double_t da,TString au,Double_t dt,TString tu,Int_t mode=1); // Provide space and time matching reference signals
73   void DisplaySignal(TString frame,TString mode,AliTimestamp* ts,Int_t jref=0,TString proj="ham",Int_t clr=0); // Display stored signal
74   void DisplaySignal(TString frame,TString mode,AliTimestamp* ts,TString name,TString proj="ham",Int_t clr=0); // Display stored signal
75   void DisplaySignals(TString frame,TString mode,AliTimestamp* ts,TString proj="ham",Int_t clr=0);             // Display all stored signals
76   void SetCentralMeridian(Double_t phi,TString u="deg"); // Set central meridian for the sky display
77  
78  protected:
79   AliPosition fLabPos;   // Position of the lab in the terrestrial longitude-latitude frame
80   Double_t fToffset;     // Lab time offset in fractional hours w.r.t. UT
81   AliSignal* fXsig;      // Signal entry for object or event studies
82   TObjArray* fRefs;      // Array holding the reference signals
83   TRotMatrix fB;         //! The frame bias matrix for conversion of ICRS to J2000 coordinates
84   Int_t fBias;           //! Initialisation flag for fB values (0=uninitialised  1=initialised)
85   TRotMatrix fP;         //! Matrix for precession correction  
86   TRotMatrix fN;         //! Matrix for nutation correction  
87   TRotMatrix fG;         //! Matrix for conversion of equatorial to galactic coordinates
88   Int_t fGal;            //! Type indicator for fG values (0=uninitialised  1=B1950  2=J2000)
89   TRotMatrix fE;         //! Matrix for conversion of equatorial to ecliptic coordinates
90   TRotMatrix fH;         //! Matrix for conversion of equatorial to horizontal coordinates
91   TRotMatrix fL;         //! Matrix for conversion of horizontal to local-frame coordinates
92   TArrayI* fIndices;     //! Storage indices of the matching reference signals
93   void SetBmatrix();                 // Set the frame bias matrix
94   void SetPmatrix(AliTimestamp* ts); // Set precession matrix for Julian date jd w.r.t. J2000.
95   void SetNmatrix(AliTimestamp* ts); // Set nutation matrix for Julian date jd w.r.t. J2000.
96   void SetGmatrix(TString mode);     // Set the equatorial to galactic conversion matrix
97   void SetEmatrix(AliTimestamp* ts); // Set the equatorial to ecliptic conversion matrix
98   void SetHmatrix(AliTimestamp* ts); // Set the equatorial to horizontal conversion matrix
99   void Precess(Ali3Vector& r,AliTimestamp* ts1,AliTimestamp* ts2); // Correct RA and decl. for earth's precession
100   void Nutate(Ali3Vector& r,AliTimestamp* ts); // Correct RA and decl. for earth's nutation
101
102   // The skymap display facilities of Garmt
103   Double_t fMeridian;  //! Central meridian (in rad) for the sky display
104   TString fProj;       //! Projection which is currently in use
105   TCanvas* fCanvas;    //! The canvas for the skymap
106   TH1* fHist;          //! Temp. histogram for the sky display
107   TObjArray* fMarkers; //! Temp. array to hold the markers for the signal display
108   void Project(Double_t l,Double_t b,TString proj,Double_t& x,Double_t& y);// Projection of (l,b) pair
109   void ProjectCylindrical(Double_t l,Double_t b,Double_t& x,Double_t& y);  // Cylindrical projection of (l,b) pair
110   void ProjectHammer(Double_t l,Double_t b,Double_t& x,Double_t& y);       // Hammer-Aitoff projection of (l,b) pair
111   void ProjectAitoff(Double_t l,Double_t b,Double_t& x,Double_t& y);       // Aitoff projection of (l,b) pair
112   void ProjectMercator(Double_t l,Double_t b,Double_t& x,Double_t& y);     // Mercator projection of (l,b) pair
113  
114  ClassDef(AliAstrolab,3) // Virtual lab to relate measurements with astrophysical phenomena
115 };
116 #endif