]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/AliMpSt345Reader.cxx
e743cd627e362b2f09811f808b130eed86e3ee80
[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.4 2005/09/19 19:01:31 ivana Exp $
18
19 #include "AliMpSt345Reader.h"
20
21 #include "AliLog.h"
22 #include "AliMpMotifReader.h"
23 #include "AliMpFiles.h"
24 #include "AliMpMotifType.h"
25 #include "AliMpPCB.h"
26 #include "AliMpSlat.h"
27 #include "AliMpMotifMap.h"
28 #include "AliMpMotifPosition.h"
29 #include "AliMpMotif.h"
30 #include "AliMpHelper.h"
31
32 #include "Riostream.h"
33
34 #include "TClass.h"
35 #include "TObjString.h"
36 #include "TString.h"
37
38 #include <sstream>
39
40 ClassImp(AliMpSt345Reader)
41
42 TMap AliMpSt345Reader::fgPCBMap;
43
44 //_____________________________________________________________________________
45 AliMpSt345Reader::AliMpSt345Reader() : TObject()
46 {
47   //
48   // Default ctor.
49   //
50
51
52 //_____________________________________________________________________________
53 AliMpSt345Reader::~AliMpSt345Reader()
54 {
55   //
56   // Dtor.
57   //
58   fgPCBMap.Delete();
59 }
60
61 //_____________________________________________________________________________
62 AliMpPCB*
63 AliMpSt345Reader::PCB(const char* pcbType)
64 {
65   //
66   // Get access to an AliMpPCB object, given its type (e.g. N1, SB2, etc...)
67   //
68   // Note that the returned object is either a new one (read from file) or a 
69   // reused one if it is already present in the internal map.
70   //
71   
72   TPair* pair = (TPair*)fgPCBMap.FindObject(pcbType);
73   if ( pair )
74   {
75     AliDebugClass(1,Form("Getting pcb %s from internal map",pcbType));
76     return (AliMpPCB*)pair->Value();
77   }
78   else
79   {
80     AliDebugClass(1,Form("Reading pcb %s from file",pcbType));
81     return ReadPCB(pcbType);
82   }
83 }
84
85 //_____________________________________________________________________________
86 AliMpPCB*
87 AliMpSt345Reader::ReadPCB(const char* pcbType)
88
89   //
90   // Create a new AliMpPCB object, by reading it from file.
91   //
92   
93   std::ifstream in(AliMpFiles::Instance()->SlatPCBFilePath(pcbType).Data());
94   if (!in.good()) 
95   {
96     AliErrorClass(Form("Cannot open file for PCB %s",pcbType));
97     return 0;
98   }
99  
100   AliMpMotifReader reader(kStation345,kNonBendingPlane); 
101   // note that the nonbending
102   // parameter is of no use for station345, as far as reading motif is 
103   // concerned, as all motifs are supposed to be in the same directory
104   // (as they are shared by bending/non-bending planes).
105      
106   char line[80];
107   
108   const TString sizeKeyword("SIZES");
109   const TString motifKeyword("MOTIF");
110   
111   AliMpPCB* pcb = 0;
112   
113   while ( in.getline(line,80) )
114   {
115     if ( line[0] == '#' ) continue;
116     
117     TString sline(line);
118     
119     if ( sline(0,sizeKeyword.Length()) == sizeKeyword )
120     {
121       std::istringstream sin(sline(sizeKeyword.Length(),
122                                    sline.Length()-sizeKeyword.Length()).Data());
123       float padSizeX = 0.0;
124       float padSizeY = 0.0;
125       float pcbSizeX = 0.0;
126       float pcbSizeY = 0.0;
127       sin >> padSizeX >> padSizeY >> pcbSizeX >> pcbSizeY;
128       assert(pcb==0);
129       pcb = new AliMpPCB(pcbType,padSizeX,padSizeY,pcbSizeX,pcbSizeY);
130     }
131     
132     if ( sline(0,motifKeyword.Length()) == motifKeyword )
133     {
134       std::istringstream sin(sline(motifKeyword.Length(),
135                                    sline.Length()-motifKeyword.Length()).Data());
136       TString sMotifType;
137       int ix;
138       int iy;
139       sin >> sMotifType >> ix >> iy;
140       
141       AliMpMotifType* motifType = 
142         reader.BuildMotifType(sMotifType.Data());
143       
144       assert(pcb!=0);
145       pcb->Add(motifType,ix,iy);
146     }
147   }
148   
149   in.close();
150   
151   fgPCBMap.Add(new TObjString(pcbType),pcb);
152   return pcb;
153 }
154
155 //_____________________________________________________________________________
156 AliMpSlat*
157 AliMpSt345Reader::ReadSlat(const char* slatType, AliMpPlaneType planeType)
158 {
159   //
160   // Create a new AliMpSlat object, by reading it from file.
161   //
162   
163   std::ifstream in(AliMpFiles::Instance()->SlatFilePath(slatType,
164                                                         planeType).Data());
165   if (!in.good()) 
166   {
167     AliErrorClass(Form("Cannot read slat from %s",
168                        AliMpFiles::Instance()->
169                        SlatFilePath(slatType,planeType).Data()));
170     return 0;
171   }
172   
173   char line[80];
174   
175   const TString pcbKeyword("PCB");
176   
177   AliMpSlat* slat = new AliMpSlat(slatType);
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,pcbKeyword.Length()) == pcbKeyword )
186     {
187       TString tmp(sline(pcbKeyword.Length()+1,sline.Length()-pcbKeyword.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 = PCB(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() != 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 slat;
216               return 0;
217       }
218       slat->Add(pcbType,manuList);
219     }
220   }
221   
222   in.close();
223   
224   return slat;
225 }
226