]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliRelAlignerKalmanArray.cxx
Updates (M. Gheata)
[u/mrichter/AliRoot.git] / PWG1 / AliRelAlignerKalmanArray.cxx
1 /**************************************************************************
2   * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3   *                                                                        *
4   * Author: The ALICE Off-line Project.                                    *
5   * Contributors are mentioned in the code where appropriate.              *
6   *                                                                        *
7   * Permission to use, copy, modify and distribute this software and its   *
8   * documentation strictly for non-commercial purposes is hereby granted   *
9   * without fee, provided that the above copyright notice appears in all   *
10   * copies and that both the copyright notice and this permission notice   *
11   * appear in the supporting documentation. The authors make no claims     *
12   * about the suitability of this software for any purpose. It is          *
13   * provided "as is" without express or implied warranty.                  *
14   **************************************************************************/
15
16 ///////////////////////////////////////////////////////////////////////////////
17 //
18 //     Data container for relative ITS-TPC alignment analysis
19 //     Holds an array of AliRelAlignerKalman objects
20 //     and takes care of merging when processing data in parallel
21 //
22 //     Origin: Mikolaj Krzewicki, Nikhef, Mikolaj.Krzewicki@cern.ch
23 //
24 //////////////////////////////////////////////////////////////////////////////
25
26 #include <TNamed.h>
27 #include <TMath.h>
28 #include <TObjArray.h>
29 #include <TCollection.h>
30 #include "AliESDEvent.h"
31 #include "AliRelAlignerKalman.h"
32 #include "AliRelAlignerKalmanArray.h"
33
34 ClassImp(AliRelAlignerKalmanArray)
35
36 //______________________________________________________________________________
37 AliRelAlignerKalmanArray::AliRelAlignerKalmanArray():
38     TNamed(),
39     fArray(new TObjArray()),
40     fSaveInterval(600), //default every 10minutes
41     fTimeMatchingTolerance(20),
42     fCurrentTimeBin(0),
43     fAligner(new AliRelAlignerKalman())
44 {
45   //ctor
46 }
47
48 //______________________________________________________________________________
49 AliRelAlignerKalmanArray::AliRelAlignerKalmanArray(const char* name):
50     TNamed(name, name),
51     fArray(new TObjArray()),
52     fSaveInterval(600), //default every 10 minutes
53     fTimeMatchingTolerance(60),
54     fCurrentTimeBin(0),
55     fAligner(new AliRelAlignerKalman())
56 {
57   //ctor
58 }
59
60 //______________________________________________________________________________
61 AliRelAlignerKalmanArray::~AliRelAlignerKalmanArray()
62 {
63   //dtor
64   fArray->SetOwner();
65   delete fArray;
66   delete fAligner;
67 }
68
69 //______________________________________________________________________________
70 Long64_t AliRelAlignerKalmanArray::Merge( TCollection* list )
71 {
72   //Merge all the arrays
73   //the merge is vertical, meaning matching entries in tree are merged
74
75   AliRelAlignerKalmanArray *arrayFromList;
76   if (!list) return 0;
77   TIter next(list);
78   while ( (arrayFromList = dynamic_cast<AliRelAlignerKalmanArray*>(next())) )
79   {
80     if (arrayFromList==this) continue;
81
82     fArray->AddAll(arrayFromList->fArray); //put all objects in one array
83
84     //do the merge
85     fArray->Sort();
86     TObjArray* tmpArray = SortedMerge(fArray);
87     tmpArray->SetOwner(kTRUE);
88
89     TObjArray* newArray = dynamic_cast<TObjArray*>(tmpArray->Clone());
90     delete fArray; //takes care of all loaded objects
91     fArray = newArray;
92
93     fArray->AddLast(arrayFromList->fAligner); //add the endofrun aligner
94   }
95
96   //TODO: this can be done better!
97   //Add own endofrun aligner and clean up
98   fArray->AddLast(fAligner);
99   fArray->Sort();
100   //TObjArray* tmpArray = SortedMerge(fArray);
101   //tmpArray->SetOwner(kTRUE);
102   //TObjArray* newArray = dynamic_cast<TObjArray*>(tmpArray->Clone());
103   //delete fArray; //takes care of all loaded objects
104   //fArray = newArray;
105
106   return fArray->GetEntriesFast();
107 }
108
109 //______________________________________________________________________________
110 TObjArray* AliRelAlignerKalmanArray::SortedMerge( TObjArray* input )
111 {
112   //Merges the adjacent aligners if close enough
113   //input needs to be already sorted
114
115   UInt_t timeStampIn;
116   AliRelAlignerKalman* alignerIn;
117   AliRelAlignerKalman* alignerOut = dynamic_cast<AliRelAlignerKalman*>(input->At(0));
118   TObjArray* output = new TObjArray();  //empty array
119   output->AddLast(alignerOut);        //first object in: copy of first input element
120
121   timeStampIn = alignerOut->GetTimeStamp();
122   SetCurrentTimeBin( timeStampIn );
123   TIter next(input);
124   while ( (alignerIn = dynamic_cast<AliRelAlignerKalman*>(next())) )
125   {
126     timeStampIn = alignerIn->GetTimeStamp();
127     if ( IsInCurrentTimeBin(timeStampIn) )
128     {
129       alignerOut->Merge(alignerIn);
130     }
131     else
132     {
133       alignerOut = alignerIn;
134       output->AddLast(alignerOut);
135     }//if
136     SetCurrentTimeBin( timeStampIn );
137   }
138   return output;
139 }
140
141 //______________________________________________________________________________
142 void AliRelAlignerKalmanArray::SetCurrentTimeBin( UInt_t timestamp )
143 {
144   //set the current timebin
145   fCurrentTimeBin = TimeBin(timestamp);
146 }
147
148 //______________________________________________________________________________
149 UInt_t AliRelAlignerKalmanArray::TimeBin( UInt_t timestamp ) const
150 {
151   return (timestamp+(fSaveInterval/2))/fSaveInterval*fSaveInterval; //it's all integers!
152 }
153
154 //______________________________________________________________________________
155 Bool_t AliRelAlignerKalmanArray::IsInCurrentTimeBin( UInt_t timestamp ) const
156 {
157   //check if timestamp is within the current timebin
158   UInt_t timeDiff = (timestamp>=fCurrentTimeBin)?timestamp-fCurrentTimeBin:
159                     fCurrentTimeBin-timestamp;
160   return (timeDiff < fTimeMatchingTolerance);
161 }
162
163 ////______________________________________________________________________________
164 //void AliRelAlignerKalmanArray::AddESDEvent( AliESDEvent* event )
165 //{
166 //  //add an AliESDEvent, take care of bookkeeping
167 //  if (!fAligner) return;
168 //  if (event->GetRunNumber() != fAligner->GetRunNumber())
169 //  {
170 //    //what to do when a new run starts
171 //  }
172 //
173 //}
174 //
175 ////______________________________________________________________________________
176 Bool_t AliRelAlignerKalmanArray::AddCosmicEvent( AliESDEvent* event )
177 {
178   if (!fAligner->AddCosmicEvent(event)) return kFALSE;
179
180   UInt_t currentTimeStamp = event->GetTimeStamp();
181   UInt_t timeFromLastBinCentre = currentTimeStamp - fCurrentTimeBin;
182   UInt_t binCentre = TimeBin(currentTimeStamp);
183   UInt_t timeFromBinCentre = currentTimeStamp-binCentre;
184   UInt_t nIntervals = timeFromLastBinCentre/fSaveInterval;
185
186   //////////////////////////////////////////////////////////////////////////////
187   //only ONE SAVE PER TIMEBIN!!!!, as close as possible to the bin center
188   //////////////////////////////////////////////////////////////////////////////
189   if ( (nIntervals == 1) &&                      //is in next time bin passed centre
190        (timeFromBinCentre < fTimeMatchingTolerance) )     //and close to it
191   {
192     AddLast(new AliRelAlignerKalman(*fAligner));
193   }
194   else if ( (nIntervals > 2) ) //TODO: don't hardwire stuff!
195   {
196     //if missed a few windows save anyway at current bin centre
197     fAligner->SetTimeStamp(binCentre);
198     AddLast(new AliRelAlignerKalman(*fAligner));
199   }
200   //////////////////////////////////////////////////////////////////////////////
201
202   return kTRUE;
203 }
204
205 //______________________________________________________________________________
206 AliRelAlignerKalmanArray::AliRelAlignerKalmanArray( const AliRelAlignerKalmanArray& in):
207     TNamed(in.GetName(), in.GetTitle()),
208     fArray(NULL),
209     fSaveInterval(in.fSaveInterval),
210     fTimeMatchingTolerance(in.fTimeMatchingTolerance),
211     fCurrentTimeBin(in.fCurrentTimeBin),
212     fAligner(new AliRelAlignerKalman(*in.fAligner))
213 {
214   //copy ctor
215   fArray = static_cast<TObjArray*>(in.Clone());
216 }
217
218 //______________________________________________________________________________
219 AliRelAlignerKalmanArray& AliRelAlignerKalmanArray::operator=(const AliRelAlignerKalmanArray& in)
220 {
221   //assignment operator
222   fArray = static_cast<TObjArray*>(in.Clone());
223   fSaveInterval = in.fSaveInterval;
224   fTimeMatchingTolerance = in.fTimeMatchingTolerance;
225   return *this;
226 }
227
228 //______________________________________________________________________________
229 Bool_t AliRelAlignerKalmanArray::SetSaveInterval( const UInt_t s )
230 {
231   //only set if array empty
232   if (fArray->GetEntriesFast()) return kFALSE;
233   fSaveInterval = s;
234   return kTRUE;
235 }
236
237 //______________________________________________________________________________
238 Bool_t AliRelAlignerKalmanArray::SetTimeMatchingTolerance( const UInt_t m )
239 {
240   //only set if array empty
241   if (fArray->GetEntriesFast()) return kFALSE;
242   fTimeMatchingTolerance = m;
243   return kTRUE;
244 }
245
246 //______________________________________________________________________________
247 AliRelAlignerKalman* AliRelAlignerKalmanArray::At( Int_t i ) const
248 {
249   //mimic TObjArray::At( Int_t i )
250   return dynamic_cast<AliRelAlignerKalman*>(fArray->At(i));
251 }
252
253 //______________________________________________________________________________
254 void AliRelAlignerKalmanArray::AddLast( AliRelAlignerKalman* al )
255 {
256   //mimic TObjArray::AddLast( TObject* obj )
257   fArray->AddLast( al );
258   SetCurrentTimeBin(al->GetTimeStamp());
259 }
260
261 //______________________________________________________________________________
262 AliRelAlignerKalman* AliRelAlignerKalmanArray::Last() const
263 {
264   //mimic TObjArray::Last()
265   return dynamic_cast<AliRelAlignerKalman*>(fArray->Last());
266 }
267
268 //______________________________________________________________________________
269 AliRelAlignerKalman* AliRelAlignerKalmanArray::operator[](Int_t i) const
270 {
271   //mimic TObjArray::operator[](Int_t)
272   return dynamic_cast<AliRelAlignerKalman*>(fArray->At(i));
273 }
274