]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/AliMpSt345Reader.cxx
bug in t0 in optimization fixed
[u/mrichter/AliRoot.git] / MUON / mapping / AliMpSt345Reader.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 purpeateose. It is      *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 // $Id$
17 // $MpId: AliMpSt345Reader.cxx,v 1.11 2006/05/24 13:58:50 ivana Exp $
18
19 #include "AliMpSt345Reader.h"
20
21 #include "AliLog.h"
22 #include "AliMpSlatMotifMap.h"
23 #include "AliMpMotifReader.h"
24 #include "AliMpFiles.h"
25 #include "AliMpMotifType.h"
26 #include "AliMpPCB.h"
27 #include "AliMpSlat.h"
28 #include "AliMpMotifPosition.h"
29 #include "AliMpMotif.h"
30 #include "AliMpHelper.h"
31 #include "AliMpConstants.h"
32
33 #include "Riostream.h"
34 #include "TClass.h"
35 #include "TObjString.h"
36 #include "TString.h"
37
38 #include <sstream>
39
40  
41 //-----------------------------------------------------------------------------
42 /// \class AliMpSt345Reader
43 //
44 /// Read slat and pcb ASCII files.
45 /// 
46 /// Basically this class provides two methods :
47 /// - AliMpSlat* ReadSlat()
48 /// - AliMpPCB ReadPCB()
49 ///
50 /// \author Laurent Aphecetche
51 //-----------------------------------------------------------------------------
52
53 /// \cond CLASSIMP
54 ClassImp(AliMpSt345Reader)
55 /// \endcond
56
57 //_____________________________________________________________________________
58 AliMpSt345Reader::AliMpSt345Reader() 
59
60 TObject(),
61 fMotifMap(AliMpSlatMotifMap::Instance())
62 {
63   ///
64   /// Default ctor.
65   ///
66
67
68 //_____________________________________________________________________________
69 AliMpSt345Reader::~AliMpSt345Reader()
70 {
71   ///
72   /// Dtor.
73   ///
74 }
75
76 //_____________________________________________________________________________
77 AliMpPCB*
78 AliMpSt345Reader::ReadPCB(const char* pcbType)
79
80   ///
81   /// Create a new AliMpPCB object, by reading it from file.
82   /// The returned object must be deleted by the client
83   
84   std::ifstream in(AliMpFiles::SlatPCBFilePath(AliMp::kStation345,pcbType).Data());
85   if (!in.good()) 
86   {
87     AliErrorClass(Form("Cannot open file for PCB %s",pcbType));
88     return 0;
89   }
90  
91   AliMpMotifReader reader(AliMp::kStation345,AliMp::kNonBendingPlane); 
92   // note that the nonbending
93   // parameter is of no use for station345, as far as reading motif is 
94   // concerned, as all motifs are supposed to be in the same directory
95   // (as they are shared by bending/non-bending planes).
96      
97   char line[80];
98   
99   const TString kSizeKeyword("SIZES");
100   const TString kMotifKeyword("MOTIF");
101   
102   AliMpPCB* pcb = 0;
103   
104   while ( in.getline(line,80) )
105   {
106     if ( line[0] == '#' ) continue;
107     
108     TString sline(line);
109     
110     if ( sline(0,kSizeKeyword.Length()) == kSizeKeyword )
111     {
112       std::istringstream sin(sline(kSizeKeyword.Length(),
113                                    sline.Length()-kSizeKeyword.Length()).Data());
114       double padSizeX = 0.0;
115       double padSizeY = 0.0;
116       double pcbSizeX = 0.0;
117       double pcbSizeY = 0.0;
118       sin >> padSizeX >> padSizeY >> pcbSizeX >> pcbSizeY;
119       if (pcb)
120       {
121         AliError("pcb not null as expected");
122       }
123       pcb = new AliMpPCB(fMotifMap,pcbType,padSizeX,padSizeY,pcbSizeX,pcbSizeY);
124     }
125     
126     if ( sline(0,kMotifKeyword.Length()) == kMotifKeyword )
127     {
128       std::istringstream sin(sline(kMotifKeyword.Length(),
129                                    sline.Length()-kMotifKeyword.Length()).Data());
130       TString sMotifType;
131       int ix;
132       int iy;
133       sin >> sMotifType >> ix >> iy;
134       
135       AliMpMotifType* motifType = fMotifMap->FindMotifType(sMotifType);
136       if (!motifType)
137       {
138         AliDebug(1,Form("Reading motifType %s from file",sMotifType.Data()));
139         motifType = reader.BuildMotifType(sMotifType.Data());
140         fMotifMap->AddMotifType(motifType);
141       }
142       else
143       {
144         AliDebug(1,Form("Got motifType %s from motifMap",sMotifType.Data()));
145       }
146       
147       pcb->Add(motifType,ix,iy);
148     }
149   }
150   
151   in.close();
152   
153   return pcb;
154 }
155
156 //_____________________________________________________________________________
157 AliMpSlat*
158 AliMpSt345Reader::ReadSlat(const char* slatType, AliMp::PlaneType planeType)
159 {
160   ///
161   /// Create a new AliMpSlat object, by reading it from file.
162   /// The returned object must be deleted by the client.
163   
164   std::ifstream in(AliMpFiles::SlatFilePath(AliMp::kStation345,slatType,
165                                             planeType).Data());
166   if (!in.good()) 
167   {
168     AliErrorClass(Form("Cannot read slat from %s",
169                        AliMpFiles::SlatFilePath(AliMp::kStation345,slatType,planeType).Data()));
170     return 0;
171   }
172   
173   char line[80];
174   
175   const TString kpcbKeyword("PCB");
176   
177   AliMpSlat* slat = new AliMpSlat(slatType, planeType);
178   
179   while ( in.getline(line,80) )
180   {
181     if ( line[0] == '#' ) continue;
182     
183     TString sline(AliMpHelper::Normalize(line));
184     
185     if ( sline(0,kpcbKeyword.Length()) == kpcbKeyword )
186     {
187       TString tmp(sline(kpcbKeyword.Length()+1,sline.Length()-kpcbKeyword.Length()));
188       Ssiz_t blankPos = tmp.First(' ');
189       if ( blankPos < 0 )
190             {
191         AliErrorClass("Syntax error in PCB file, should get a list of "
192                       "manu ids after the pcbname");
193         delete slat;
194         return 0;
195             }
196       
197       TString pcbName(tmp(0,blankPos));
198       TString manus(tmp(blankPos+1,tmp.Length()-blankPos));
199       
200       AliMpPCB* pcbType = ReadPCB(pcbName.Data());        
201       if (!pcbType)
202             {
203         AliErrorClass(Form("Cannot read pcbType=%s",pcbName.Data()));
204               delete slat;
205               return 0;
206             }      
207
208       TArrayI manuList;
209       AliMpHelper::DecodeName(manus,';',manuList);
210       if ( manuList.GetSize() != Int_t(pcbType->GetSize()) )
211             {
212         AliErrorClass(Form("Wrong number of manu ids for this PCB ("
213                            "%s) : %d out of %d",pcbName.Data(),
214                            manuList.GetSize(),pcbType->GetSize()));
215         delete pcbType;
216               delete slat;
217               return 0;
218       }
219
220       for ( Int_t i = 0; i < manuList.GetSize(); ++i )
221       {
222         manuList[i] |= AliMpConstants::ManuMask(planeType);
223       }
224       slat->Add(*pcbType,manuList);
225       delete pcbType;
226     }
227   }
228   
229   in.close();
230   
231   return slat;
232 }
233