]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPedestalSubprocessor.cxx
Includes required by ROOT head
[u/mrichter/AliRoot.git] / MUON / AliMUONPedestalSubprocessor.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 // $Id$
17
18 #include <sstream>
19
20 #include <Riostream.h>
21 #include <TList.h>
22 #include <TObjString.h>
23 #include <TSystem.h>
24
25 #include "AliCDBMetaData.h"
26 #include "AliLog.h"
27 #include "AliMUON2DMap.h"
28 #include "AliMUON2DStoreValidator.h"
29 #include "AliMUONCalibParam2F.h"
30 #include "AliMUONConstants.h"
31 #include "AliMUONObjectPair.h"
32 #include "AliMUONPedestalSubprocessor.h"
33 #include "AliMUONPreprocessor.h"
34 #include "AliMUONVDataIterator.h"
35 #include "AliMpDDLStore.h"
36
37 ///
38 /// \class AliMUONPedestalSubprocessor
39 ///
40 /// Implementation of AliMUONVSubprocessor class to deal with MUON TRK pedestals.
41 ///
42 /// Pedestals are read in from an ascii file, with the format :               \n
43 ///---------------------------------------------------------------------------\n
44 /// BUS_PATCH MANU_ADDR CHANNEL      MEAN       SIGMA                         \n
45 ///---------------------------------------------------------------------------\n
46 ///
47 /// \author L. Aphecetche
48 ///
49
50 /// \cond CLASSIMP
51 ClassImp(AliMUONPedestalSubprocessor)
52 /// \endcond
53
54 //_____________________________________________________________________________
55 AliMUONPedestalSubprocessor::AliMUONPedestalSubprocessor(AliMUONPreprocessor* master)
56 : AliMUONVSubprocessor(master,
57                        "Pedestals",
58                        "Upload MUON Tracker pedestals to OCDB"),
59 fPedestals(0x0)
60 {
61   // default ctor
62 }
63
64 //_____________________________________________________________________________
65 AliMUONPedestalSubprocessor::~AliMUONPedestalSubprocessor()
66 {
67   // dtor
68   delete fPedestals;
69 }
70
71 //_____________________________________________________________________________
72 void 
73 AliMUONPedestalSubprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
74 {
75   /// When starting a new run, reads in the pedestals ASCII files.
76   
77   const Int_t kSystem = AliMUONPreprocessor::kDAQ;
78   const char* kId = "PEDESTALS";
79   
80   delete fPedestals;
81   fPedestals = new AliMUON2DMap(kTRUE);
82   
83   AliInfo(Form("Reading pedestal files for Run %d startTime %ld endTime %ld",
84                run,startTime,endTime));
85   
86   TList* sources = Master()->GetFileSources(kSystem,kId);
87   TIter next(sources);
88   TObjString* o(0x0);
89   while ( ( o = static_cast<TObjString*>(next()) ) )
90   {
91     TString fileName(Master()->GetFile(kSystem,kId,o->GetName()));
92     Bool_t ok = ReadFile(fileName.Data());
93     if (!ok)
94     {
95       AliError(Form("Could not read file %s",fileName.Data()));
96     }
97   }
98   delete sources;
99 }
100
101 //_____________________________________________________________________________
102 UInt_t 
103 AliMUONPedestalSubprocessor::Process(TMap* /*dcsAliasMap*/)
104 {
105   // Store the pedestals into the CDB
106   
107   if (!fPedestals) return 0;
108   
109   AliInfo("Validating pedestals");
110   AliMUON2DStoreValidator validator;
111   TObjArray* missing =
112     validator.Validate(*fPedestals,AliMUONCalibParam2F::InvalidFloatValue());  
113   
114   if (missing)  
115   {
116     validator.Report(*missing);
117 //    AliError("Will not write into CDB as some pieces are missing...");
118 //    return 0;
119   }
120   
121   AliInfo("Storing pedestals");
122   
123   AliCDBMetaData metaData;
124         metaData.SetBeamPeriod(0);
125         metaData.SetResponsible("MUON TRK");
126         metaData.SetComment("Computed by AliMUONPedestalSubprocessor $Id$");
127   
128         UInt_t result = Master()->Store("Calib", "Pedestals", fPedestals, &metaData, 0, 0);
129   
130   return result;  
131 }
132
133 //_____________________________________________________________________________
134 Bool_t
135 AliMUONPedestalSubprocessor::ReadFile(const char* filename)
136 {
137   // Read the pedestals from an ASCII file.
138   // Format of that file is one line per channel :
139   //---------------------------------------------------------------------------
140   // BUS_PATCH MANU_ADDR CHANNEL      MEAN       SIGMA
141   //---------------------------------------------------------------------------
142   //
143   // Return kFALSE if reading was not successfull.
144   //
145   
146   AliInfo(Form("Reading %s",filename));
147   
148   std::ifstream in(gSystem->ExpandPathName(filename));
149   if (!in.good()) return kFALSE;
150   
151   char line[80];
152   Int_t busPatchID, manuID, manuChannel;
153   Float_t pedMean, pedSigma;
154   static const Int_t kNchannels(64);
155   static Bool_t replace(kFALSE);
156   
157   while ( in.getline(line,80) )
158   {
159     if ( line[0] == '/' && line[1] == '/' ) continue;
160     std::istringstream sin(line);
161     sin >> busPatchID >> manuID >> manuChannel >> pedMean >> pedSigma;
162     Int_t detElemID = AliMpDDLStore::Instance()->GetDEfromBus(busPatchID);
163     AliDebug(3,Form("BUSPATCH %3d DETELEMID %4d MANU %3d CH %3d MEAN %7.2f SIGMA %7.2f",
164              busPatchID,detElemID,manuID,manuChannel,pedMean,pedSigma));
165     
166     AliMUONVCalibParam* ped = 
167       static_cast<AliMUONVCalibParam*>(fPedestals->Get(detElemID,manuID));
168     
169     if (!ped) 
170     {
171       ped = new AliMUONCalibParam2F(kNchannels,AliMUONCalibParam2F::InvalidFloatValue());  
172       fPedestals->Set(detElemID,manuID,ped,replace);
173     }
174     ped->SetValueAsFloat(manuChannel,0,pedMean);
175     ped->SetValueAsFloat(manuChannel,1,pedSigma);
176   }
177   in.close();
178   return kTRUE;
179 }
180
181
182 //_____________________________________________________________________________
183 void
184 AliMUONPedestalSubprocessor::Print(Option_t* opt) const
185 {
186   /// ouput to screen
187   AliMUONVDataIterator* it = fPedestals->Iterator();
188   AliMUONObjectPair* p;
189
190   while ( ( p = static_cast<AliMUONObjectPair*>(it->Next() ) ) )
191   {
192     AliMUONVCalibParam* value = static_cast<AliMUONVCalibParam*>(p->Value());
193     value->Print(opt);
194   }
195 }