]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Adding the AliTPCkalmanTime class (Marian)
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 Dec 2008 13:39:52 +0000 (13:39 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 Dec 2008 13:39:52 +0000 (13:39 +0000)
TPC/AliTPCkalmanTime.cxx [new file with mode: 0644]
TPC/AliTPCkalmanTime.h [new file with mode: 0644]
TPC/TPCbaseLinkDef.h
TPC/libTPCbase.pkg

diff --git a/TPC/AliTPCkalmanTime.cxx b/TPC/AliTPCkalmanTime.cxx
new file mode 100644 (file)
index 0000000..366e58d
--- /dev/null
@@ -0,0 +1,143 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/*
+  TPC Kalman filter implementation for time dependent variables:
+
+
+The drift velocity and the gas gain are changing in time.
+The drift velocity and gas gain is a function of many parameters, but not all of 
+them are known. We assume that the most important parameters are pressure and temperature
+and the influence of other parameters (gas composition, and electric field) are only 
+slowly varying in time and can be expressed by smooth function $x_{off}(t)$:
+\begin{equation}
+x(t) = x_{off}(t)+k_N\frac{\Delta{P/T}}{P/T}   
+\end{equation}
+where x(t) is the parameter which we observe.
+\begin{equation}
+\begin{split}
+x(t)=\frac{\Delta{G}}{G_0}     \\
+x(t)=\frac{\Delta{v_d}}{v_{d0}}        
+\end{split}
+\end{equation}
+
+Kalman filter parameters are following:
+\begin{itemize}
+\item State vector  ($x_{off}(t)$, $k_N$) at given time
+\item Covariance matrix
+\end{itemize}
+
+Kalman filter implent following functions:
+\begin{itemize}
+\item Prediction - adding covariance element $\sigma_{xoff}$
+\item Update state vector with new measurement vector ($x_t,\frac{\Delta{P/T}}{P/T}$)
+\end{itemize}
+
+
+*/
+
+#include "AliTPCkalmanTime.h"
+
+
+AliTPCkalmanTime::AliTPCkalmanTime():
+  TNamed(),
+  fState(0),
+  fCovariance(0),
+  fTime(0)
+{
+  //
+  // Default constructor
+  //
+}
+
+AliTPCkalmanTime::AliTPCkalmanTime(Double_t time, Double_t xoff, Double_t k, Double_t sigmaxoff, Double_t sigmak): 
+  TNamed(),
+  fState(0),
+  fCovariance(0),
+  fTime(0)
+{
+  //
+  // Default constructor
+  //
+  Init(time,xoff,k,sigmaxoff,sigmak);
+}
+
+
+void AliTPCkalmanTime::Init(Double_t time, Double_t xoff, Double_t k, Double_t sigmaxoff, Double_t sigmak){
+  //
+  // Default constructor
+  //
+  fState = new TMatrixD(3,1);
+  fCovariance = new TMatrixD(2,2);
+  (*fState)(0,0)= xoff;  // offset
+  (*fState)(1,0)= k;     // slope of the taylor
+  fTime=time;
+  (*fCovariance)(0,0)=sigmaxoff*sigmaxoff;
+  (*fCovariance)(1,1)=sigmak*sigmak;
+  (*fCovariance)(0,1)=0;
+  (*fCovariance)(1,0)=0;
+}
+
+
+void AliTPCkalmanTime::Propagate(Double_t time, Double_t sigma){
+  //
+  // Propagate the Kalman
+  //
+  if (!fCovariance) return;
+  if (!fState) return;
+  Double_t deltaT  =time-(*fState)(2,0);  //delta time - param2 is the current time
+  Double_t sigmaT2 =(deltaT*deltaT)*sigma*sigma; 
+  (*fCovariance)(0,0)+=sigmaT2;
+  fTime=time;  
+}
+
+void  AliTPCkalmanTime::Update(Double_t x, Double_t xerr, Double_t ptratio){
+  //
+  //
+  //
+  static TMatrixD matHk(1,2);    // vector to mesurement
+  static TMatrixD measR(1,1);    // measurement error 
+  static TMatrixD matQk(2,2);    // prediction noise vector
+  static TMatrixD vecZk(1,1);    // measurement
+  //
+  static TMatrixD vecYk(1,1);    // Innovation or measurement residual
+  static TMatrixD matHkT(2,1);
+  static TMatrixD matSk(1,1);    // Innovation (or residual) covariance
+  static TMatrixD matKk(2,1);    // Optimal Kalman gain
+  static TMatrixD mat1(2,2);     // update covariance matrix
+  static TMatrixD covXk2(2,2);
+  //
+  //
+  TMatrixD covXk(*fCovariance);    // X covariance 
+  TMatrixD vecXk(*fState);         // X vector
+  //
+  vecZk(0,0) = x;                // current mesurement
+  measR(0,0) = xerr*xerr;        // measurement error                      
+  matHk(0,0)=1;  matHk(0,1)= ptratio; 
+  vecYk = vecZk-matHk*vecXk;     // Innovation or measurement residual
+  //
+  //
+  matHkT=matHk.T(); matHk.T();
+  matSk = (matHk*(covXk*matHkT))+measR;    // Innovation (or residual) covariance
+  matSk.Invert();
+  matKk = (covXk*matHkT)*matSk;            //  Optimal Kalman gain
+  vecXk += matKk*vecYk;                    //  updated vector 
+  mat1(0,0)=1; mat1(1,1)=1; mat1(2,2)=1;
+  covXk2= (mat1-(matKk*matHk));
+  covXk =  covXk2*covXk;  
+  //
+  (*fCovariance)=covXk;
+  (*fState)=vecXk;
+}
diff --git a/TPC/AliTPCkalmanTime.h b/TPC/AliTPCkalmanTime.h
new file mode 100644 (file)
index 0000000..56bd43f
--- /dev/null
@@ -0,0 +1,26 @@
+#ifndef ALITPCKALMANTIME_H
+#define ALITPCKALMANTIME_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+#include "TNamed.h"
+#include "TMatrixD.h"
+
+
+class AliTPCkalmanTime: public TNamed{
+public:
+  AliTPCkalmanTime();
+  AliTPCkalmanTime(Double_t time, Double_t xoff, Double_t k, Double_t sigmaxoff, Double_t sigmak);
+  void Propagate(Double_t time, Double_t sigma);
+  void Update(Double_t x, Double_t xerr, Double_t ptratio);
+public:
+  void Init(Double_t time, Double_t xoff, Double_t k, Double_t sigmaxoff, Double_t sigmak);
+  TMatrixD * fState;           // state vector
+  TMatrixD * fCovariance;      // covariance
+  Double_t   fTime;            // current time
+  ClassDef(AliTPCkalmanTime,1);
+};
+
+#endif
+
index 3587cc1b0ff373ae44f0b9f3753a7e9a31cab390..f9c8c823a52c9eb7bd38e11c1791e1131121f062 100644 (file)
@@ -73,5 +73,6 @@
 #pragma link C++ class AliTPCdataQA;
 #pragma link C++ class AliTPCQAChecker+;
 #pragma link C++ class AliExternalComparison+;
+#pragma link C++ class AliTPCkalmanTime+;
 #endif
 
index 56516926e48f5d39ad2ec8361149a733e30d6c18..cfea08afebe0118fe3d73a590752a5cda9895791 100644 (file)
@@ -19,7 +19,8 @@ SRCS:=  AliSegmentID.cxx  AliSegmentArray.cxx AliDigits.cxx AliH2F.cxx \
         AliTPCTempMap.cxx AliTPCCalibVdrift.cxx \
        AliTPCCalibTCF.cxx AliTPCAltroEmulator.cxx  AliTPCeventInfo.cxx \
         AliTransform.cxx AliTPCTransform.cxx AliTPCAlign.cxx AliTPCGoofieValues.cxx \
-       AliTPCdataQA.cxx AliTPCQAChecker.cxx AliTPCConfigDA.cxx AliExternalComparison.cxx
+       AliTPCdataQA.cxx AliTPCQAChecker.cxx AliTPCConfigDA.cxx AliExternalComparison.cxx \
+       AliTPCkalmanTime.cxx
 
 HDRS:= $(SRCS:.cxx=.h)
 HDRS:= $(filter-out AliTPCTrackHitsInterfaces.h, $(HDRS))