]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDRawStream.h
Call to HTA added.
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDRawStream.h
1 #ifndef ALIHMPIDRAWSTREAM_H
2 #define ALIHMPIDRAWSTREAM_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 ///////////////////////////////////////////////////////////////////////////////
7 ///
8 /// This is a class for reading raw data digits for HMPID.
9 /// The data format is taken from the document provided by Paolo Martinengo.
10 ///
11 /// cvetan.cheshkov@cern.ch 19/07/2007
12 ///
13 ///////////////////////////////////////////////////////////////////////////////
14
15 #include <TObject.h>
16 #include <TRandom.h>
17 #include "AliHMPIDParam.h"
18 #include <AliBitPacking.h>
19 #include <AliFstream.h>
20 #include "AliHMPIDDigit.h"
21 #include "AliDAQ.h"
22 class AliRawReader;
23
24 class AliHMPIDRawStream: public TObject {
25   public :
26     AliHMPIDRawStream(AliRawReader* rawReader);
27     AliHMPIDRawStream();
28     
29     virtual ~AliHMPIDRawStream();
30
31     virtual void     Reset();
32     virtual Bool_t   Next();
33             void     Init();
34     
35             Int_t    Ch( Int_t ddl,Int_t row,Int_t dil,Int_t pad    ) {return AliHMPIDParam::A2C(fPad[ddl][row][dil][pad]); }            //chamber number
36             
37             Int_t GetDDLNumber()  const { return fDDLNumber; }                                                                            // Provide current DDL number
38     inline Int_t GetCharge(Int_t ddl,Int_t row, Int_t dilogic, Int_t pad);                                                                         // Provide the charge observed in certain row,dilogic,pad channel
39             inline Int_t GetPad(Int_t ddl,Int_t row,Int_t dil,Int_t pad);                                                                        //
40             
41             Int_t   Pc          ( Int_t ddl,Int_t row,Int_t dil,Int_t pad                            ) {return AliHMPIDParam::A2P(fPad[ddl][row][dil][pad]);}                                                 //PC position number
42             Int_t   PadPcX      ( Int_t ddl,Int_t row,Int_t dil,Int_t pad                            ) {return AliHMPIDParam::A2X(fPad[ddl][row][dil][pad]);}                                                 //pad pc x # 0..79
43             Int_t   PadPcY      ( Int_t ddl,Int_t row,Int_t dil,Int_t pad                            ) {return AliHMPIDParam::A2Y(fPad[ddl][row][dil][pad]);}                                                 //pad pc y # 0..47
44    
45             inline  Bool_t SetZeroSup (Bool_t isSup);
46     inline  Bool_t GetZeroSup(); 
47                  
48     inline void    Raw            (UInt_t &w32,Int_t &ddl,Int_t &r,Int_t &d,Int_t &a);                                               //digit->(w32,ddl,r,d,a)
49     inline void    Raw            (Int_t ddl,Int_t r,Int_t d,Int_t a);                                                               //raw->abs pad number
50     inline Bool_t  Raw            (UInt_t  w32,Int_t  ddl,AliRawReader *pRR);                                                      //(w32,ddl)->digit
51     inline void   SetCharge      (Int_t ddl,Int_t row,Int_t dil,Int_t pad,Int_t q);
52     inline void    WriteRaw       (TObjArray *pDigLst                             );                                                      //write as raw stream     
53     inline void   WriteRowMarker  (AliFstream *ddl,UInt_t size);
54     inline void   WriteEoE        (AliFstream *ddl,UInt_t row,UInt_t dil,UInt_t wordCnt);  
55     inline void   WriteSegMarker  (AliFstream *ddl,UInt_t row);   
56     
57 //    inline TClonesArray  ReMap(TClonesArray *pDigIn);
58     
59     enum EHMPIDRawStreamError {
60       kRawDataSizeErr = 1,
61       kRowMarkerErr = 2,
62       kWrongRowErr = 3,
63       kWrongDilogicErr = 4,
64       kWrongPadErr = 5,
65       kEoEFlagErr = 6,
66       kEoESizeErr = 7,
67       kEoEDILOGICErr = 8,
68       kEoERowErr = 9,
69       kBadSegWordErr = 10,
70       kWrongSegErr = 11
71     };
72     
73     enum {
74       kNRows       = 24,                                    // Number of rows (starting from 1 !)//was25
75       kNDILOGICAdd = 10,                                    // Number of DILOGIC addresses in a row (starting from 1 !) //was11
76       kNPadAdd     = 48,                                    // Number of pad row
77       kNRowsPerSegment = 8,                                 // Number of rows per segment
78       kNDDL = 14
79     };
80    enum EHMPIDRawError {
81     kInvalidRawDataWord = 1
82   };
83
84     
85   private :
86
87     AliHMPIDRawStream& operator = (const AliHMPIDRawStream& stream);
88     AliHMPIDRawStream(const AliHMPIDRawStream& stream);
89
90     UInt_t           GetNextWord();
91
92     Int_t            fCharge[kNDDL][kNRows+1][kNDILOGICAdd+1][kNPadAdd]; // Array for charge values for all channels in one DDL
93     
94     Int_t            fPad[kNDDL][kNRows+1][kNDILOGICAdd+1][kNPadAdd]; // Array for abs pad values for all channels in one DDL
95
96     UInt_t           fRawWord[kNDDL][kNRows+1][kNDILOGICAdd+1][kNPadAdd];
97         
98     Int_t            fDDLNumber;    // index of current DDL number
99
100     AliRawReader*    fRawReader;    // object for reading the raw data
101
102     UChar_t*         fData;         // raw data
103
104     Int_t            fPosition;     // current position in fData
105     
106     Bool_t           fZeroSup;
107
108     ClassDef(AliHMPIDRawStream, 0)  // base class for reading HMPID raw digits
109 };
110 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111 void AliHMPIDRawStream::Raw(UInt_t &w32,Int_t &ddl,Int_t &r,Int_t &d,Int_t &a)
112 {
113 // Convert raw stream word to raw word format
114 // Arguments: w32,ddl,r,d,a where to write the results
115 // Returns: none
116   Int_t y2a[6]={5,3,1,0,2,4};
117   
118   ddl=2*Ch(ddl,r,d,a)+Pc(ddl,r,d,a)%2;                                                          //DDL# 0..13
119   Int_t tmp=1+Pc(ddl,r,d,a)/2*8+PadPcY(ddl,r,d,a)/6;  r=(Pc(ddl,r,d,a)%2)? 25-tmp:tmp;              //row r=1..24
120   d=1+PadPcX(ddl,r,d,a)/8;                                                                  //DILOGIC# 1..10
121   a=y2a[PadPcY(ddl,r,d,a)%6]+6*(PadPcX(ddl,r,d,a)%8);                                           //ADDRESS 0..47        
122   
123   w32=0;    
124   AliBitPacking::PackWord((fCharge[ddl][r][d][a]>4095)?4095:(UInt_t)fCharge[ddl][r][d][a],w32, 0,11);       // 0000 0rrr rrdd ddaa aaaa qqqq qqqq qqqq        Qdc               bits (00..11) counts (0..4095)
125   //molnarl: Since in simulation the the charge can be > than 4095 but not in real life we need to protect. If fQ>4095 after packing we will get 0 for the charge! 
126   assert(0<=a&&a<=47);AliBitPacking::PackWord(        a ,w32,12,17);  // 3322 2222 2222 1111 1111 1000 0000 0000        DILOGIC address   bits (12..17) counts (0..47)
127   assert(1<=d&&d<=10);AliBitPacking::PackWord(        d ,w32,18,21);  // 1098 7654 3210 9876 5432 1098 7654 3210        DILOGIC number    bits (18..21) counts (1..10)
128   assert(1<=r&&r<=24);AliBitPacking::PackWord(        r ,w32,22,26);  //                                                Row number        bits (22..26) counts (1..24)  
129 }
130 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
131 void AliHMPIDRawStream::Raw(Int_t ddl,Int_t r,Int_t d,Int_t a)
132 {
133   //Assign absolute pad ID based on ddl,row,dil,pad
134   //Arguments: DDL, row number, dilogic number, dilogic address(pad)
135   //Returns  : nothing
136
137   assert(0<=ddl&&ddl<=13); assert(1<=r&&r<=24); assert(1<=d&&d<=10);   assert(0<=a&&a<=47);  
138   Int_t a2y[6]={3,2,4,1,5,0};//pady for a given address (for single DILOGIC chip)
139                                   Int_t ch=ddl/2;
140   Int_t tmp=(r-1)/8;              Int_t pc=(ddl%2)? 5-2*tmp:2*tmp; 
141                                   Int_t px=(d-1)*8+a/6;
142         tmp=(ddl%2)?(24-r):r-1;   Int_t py=6*(tmp%8)+a2y[a%6];
143   fPad[ddl][r][d][a]=AliHMPIDParam::Abs(ch,pc,px,py);
144 }
145 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
146 Bool_t AliHMPIDRawStream::Raw(UInt_t w32,Int_t ddl, AliRawReader *pRR)
147 {
148 // Converts a given raw data word to a digit
149 // Arguments: w32 - 32 bits raw data word
150 //            ddl - DDL idx  0 1 2 3 4 ... 13
151 //   Returns: none
152   Int_t r = AliBitPacking::UnpackWord(w32,22,26); assert(1<=r&&r<=24);   //                                         Row number      (1..24)    
153   Int_t d = AliBitPacking::UnpackWord(w32,18,21); assert(1<=d&&d<=10);   // 3322 2222 2222 1111 1111 1000 0000 0000 DILOGIC number  (1..10)
154   Int_t a = AliBitPacking::UnpackWord(w32,12,17); assert(0<=a&&a<=47);   // 1098 7654 3210 9876 5432 1098 7654 3210 DILOGIC address (0..47)  
155   Int_t q = AliBitPacking::UnpackWord(w32, 0,11); assert(0<=q&&q<=4095); // 0000 0rrr rrdd ddaa aaaa qqqq qqqq qqqq Qdc             (0..4095) 
156   if (r<1 || r>24 || d<1 || d>10 || a<0 || a>47 || q<0 || q>4095) {
157     AliWarning(Form("Invalid raw data word %x",w32));
158     pRR->AddMajorErrorLog(kInvalidRawDataWord,Form("w=%x",w32));
159     return kFALSE;
160   }
161   Raw(ddl,r,d,a);
162   fCharge[ddl][r][d][a]=q;
163   return kTRUE;
164 }
165 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
166 void AliHMPIDRawStream::SetCharge(Int_t ddl,Int_t row,Int_t dil,Int_t pad,Int_t q)
167
168   //Setter for the charger in the raw stream
169   //Arguments: DDL, row number, dilogic number, dilogic address(pad), charge
170   //Returns:   Charge from the raw stream
171      fCharge[ddl][row][dil][pad]=q;
172      //  return fCharge[ddl][row][dil][pad];
173       
174 }
175 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
176 Int_t AliHMPIDRawStream::GetPad(Int_t ddl,Int_t row,Int_t dil,Int_t pad)
177 {
178   // The method returns the absolute pad number or -1 
179   // in case the charge from the channels
180   // has not been read or invalid arguments
181  
182   assert(0<=ddl&&ddl<=13);  
183   assert(1<=row&&row<=24); 
184   assert(1<=dil&&dil<=10);   
185   assert(0<=pad&&pad<=47);  
186   
187   Int_t a2y[6]={3,2,4,1,5,0};     //pady for a given padress (for single DILOGIC chip)
188                                   Int_t ch=ddl/2;
189   Int_t tmp=(row-1)/8;            Int_t pc=(ddl%2)? 5-2*tmp:2*tmp; 
190                                   Int_t px=(dil-1)*8+pad/6;
191         tmp=(ddl%2)?(24-row):row-1;   
192                                   Int_t py=6*(tmp%8)+a2y[pad%6];
193                                   
194   return fPad[ddl][row][dil][pad]=AliHMPIDParam::Abs(ch,pc,px,py);
195 }//GetPad()
196 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
197 Int_t AliHMPIDRawStream::GetCharge(Int_t ddl,Int_t row, Int_t dilogic, Int_t pad)
198 {
199   // The method returns the charge collected
200   // in a particular channel
201   // Return -1 in case the charge from the channels
202   // has not been read or invalid arguments
203   if (ddl < 0 || ddl > kNDDL) {
204     AliError(Form("Wrong DDL index %d!",ddl));
205     return 0;
206   }  
207   if (row < 1 || row > kNRows) {
208     AliError(Form("Wrong row index %d!",row));
209     return 0;
210   }
211
212   if (dilogic < 1 || dilogic > kNDILOGICAdd) {
213     AliError(Form("Wrong DILOGIC address %d!",dilogic));
214     return 0;
215   }
216   
217   if (pad >= kNPadAdd) {
218     AliError(Form("Wrong pad index %d!",pad));
219     return 0;
220   }
221
222   return fCharge[ddl][row][dilogic][pad];
223   
224 }
225 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
226 void AliHMPIDRawStream::WriteRowMarker(AliFstream *ddl,UInt_t size)
227 {
228   //Writes the row marker for real data and pedestal into the ddl stream
229   //Arguments: ddl stream and the size of the block of the given row, the siye is at least the 10 EoE words!
230   //Returns:   nothing
231   UInt_t w32=0;
232   UInt_t marker=12968;                                    //ror marker: 32a8 in hexa; 12968 in decimal
233   AliBitPacking::PackWord(size,  w32, 16,31);            //number of roaw written after row marker (digits and EoE)
234   AliBitPacking::PackWord(marker,w32,0,15);              //32a8=12968
235   ddl->WriteBuffer((char*)&w32,sizeof(w32));              
236 }
237 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
238 void AliHMPIDRawStream::WriteEoE(AliFstream *ddl,UInt_t row,UInt_t dil,UInt_t wordCnt  )
239 {
240   //Writes the EoE word from real data and pedestals into the ddl stream
241   //Arguments:  ddl stream, row number, dilogic number and the number of words before the EoE
242   //Retursns:   nothing
243   UInt_t e=1;
244   UInt_t w32=0;
245   assert(1<=row&&row<=24);      AliBitPacking::PackWord((UInt_t)row     ,w32,22,26);    // row number (1...24)
246   assert(1<=dil&&dil<=10);      AliBitPacking::PackWord((UInt_t)dil     ,w32,18,21);    // DILOGIC number (1...10)
247                                 AliBitPacking::PackWord(          e     ,w32, 7,17);   // event number -- not used
248                                 AliBitPacking::PackWord((UInt_t)wordCnt ,w32, 0, 6);  // word counter (0...47)                                                                  AliBitPacking::PackWord((UInt_t)1       ,w32,27,27);  // bit 27 is always 1 by definition of EoE
249                                 AliBitPacking::PackWord((UInt_t)1       ,w32,27,27);  // bit 27 is always 1 by definition of EoE    
250   ddl->WriteBuffer((char*)&w32,sizeof(w32));      
251
252 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++     
253 void AliHMPIDRawStream::WriteSegMarker(AliFstream *ddl,UInt_t row)
254 {
255   //Writes the segment marker (after 8 rows) into the ddl stream
256   //Arguments: ddl stream and the segment: row 8 -> 0x5800, row 16 -> 5801, row 24 -> 5802
257   //Retruns:   nothing
258     UInt_t w32=0;
259                                 AliBitPacking::PackWord((UInt_t)43791,        w32,16,31);       //43791==AB0F
260                                 AliBitPacking::PackWord((UInt_t)(22528+row/8),w32, 0,15);       //22528==5800
261     
262     ddl->WriteBuffer((char*)&w32,sizeof(w32)); 
263 }      
264 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++     
265 Bool_t AliHMPIDRawStream::SetZeroSup (Bool_t isSup)
266 {
267   //Prevision to turn OFF zero suppression
268   //Arguments: setter
269   //Returns:   switch
270   fZeroSup=isSup;
271   return fZeroSup;
272 }
273 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++     
274 Bool_t AliHMPIDRawStream::GetZeroSup()
275 {
276   if(fZeroSup==kTRUE) return kTRUE;
277   else                return kFALSE;
278 }
279 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++     
280 void AliHMPIDRawStream::WriteRaw(TObjArray *pDigAll)
281 {
282 // Write a list of digits for a given chamber in raw data stream
283 // Arguments: pDigAll- list of digits 
284 //   Returns: none      
285   Int_t  ddl,r,d,a;            //32 bits data word 
286   Int_t  cntLpad,cntRpad;
287   Int_t  cntLrow,cntRrow;
288   Int_t  cntL=0,cntR=0;                           //data words counters for DDLs
289   Int_t  cntLeoe,cntReoe;
290   UInt_t posL,posR;
291   UInt_t cntLseg,cntRseg;
292   Int_t  cntRdig=0,cntLdig=0;
293   
294   UInt_t posLmarker,posRmarker;
295   Int_t digcnt=0;
296
297   Int_t isDigThere[14][25][11][48];
298   
299   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//chambers loop
300     cntL=0;cntR=0;   
301     for(Int_t iddl=0;iddl<14;iddl++){
302       for(Int_t irow=1;irow<=24;irow++){
303         for(Int_t idil=1;idil<=10;idil++){
304           for(Int_t ipad=0;ipad<48;ipad++){
305             isDigThere[iddl][irow][idil][ipad]=-1;
306           }
307         }
308       }
309     }
310     
311     AliFstream* ddlL;                                 //output streams, 2 per chamber
312     AliFstream* ddlR;                          
313     
314     AliRawDataHeader header; header.SetAttribute(0);  //empty DDL header
315     
316     ddlL = new AliFstream(AliDAQ::DdlFileName("HMPID",2*iCh+1)); //left and right looking at the IP
317     ddlR = new AliFstream(AliDAQ::DdlFileName("HMPID",2*iCh));   //open both DDL of this chamber in parallel
318     
319     ddlL->WriteBuffer((char*)&header,sizeof(header));            //write dummy header as place holder, actual 
320     ddlR->WriteBuffer((char*)&header,sizeof(header));            //will be rewritten later when total size of DDL is known
321     
322     UInt_t w32=0;                 //32 bits data word 
323     digcnt=0;
324     
325     TClonesArray *pDigCh=(TClonesArray *)pDigAll->At(iCh); //list of digits for current chamber 
326     //Printf("::::::::::::::::::; Number of Digits to write: %d == %d",pDigCh->GetEntriesFast(),pDigCh->GetEntries());
327     for(Int_t iDig=0;iDig<pDigCh->GetEntriesFast();iDig++){//digits loop
328       AliHMPIDDigit *pDig1=(AliHMPIDDigit*)pDigCh->At(iDig);
329       pDig1->Raw(w32,ddl,r,d,a);
330       isDigThere[ddl][r][d][a]=iDig;
331     }  
332     
333     for(Int_t row = 1; row <= AliHMPIDRawStream::kNRows; row++){ //AliHMPIDRawStream::kNRows=25!
334       cntRrow=0;cntLrow=0;cntLseg=0;cntRseg=0;// 
335       cntLeoe=0;cntReoe=0;
336       posLmarker=ddlL->Tellp(); WriteRowMarker(ddlL,(UInt_t)1);   cntL++; cntRrow++; 
337       posRmarker=ddlR->Tellp(); WriteRowMarker(ddlR,(UInt_t)1);   cntR++; cntLrow++; 
338       for(Int_t dil = 1; dil <= AliHMPIDRawStream::kNDILOGICAdd; dil++){ //AliHMPIDRawStream::kNDILOGICAdd = 11!
339         cntLpad=0;cntRpad=0;
340         for(Int_t pad = 0; pad < AliHMPIDRawStream::kNPadAdd; pad++){   //AliHMPIDRawStream::kNPadAdd     = 48
341           for ( Int_t iddl=2*iCh; iddl<=2*iCh+1;iddl++){
342             if (isDigThere[iddl][row][dil][pad]!=-1) {
343               AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCh->At(isDigThere[iddl][row][dil][pad]);             
344               pDig->Raw(w32,ddl,r,d,a);  
345               if(pDig->Q() < 0 ) continue;                                                 //We can turn of the zero sup for pedestal simulation
346               //Printf("::::::::::::::: ddl from Digit : %d",ddl);
347               if(ddl%2){                                                                               //write raw digit selecting on DDL
348                 ddlL->WriteBuffer((char*)&w32,sizeof(w32));   cntL++; cntLpad++; cntLrow++;  cntLdig++;//Printf(" WL: %x isDig: %d",w32,isDigThere[iddl][row][dil][pad]);
349               }else{
350                 ddlR->WriteBuffer((char*)&w32,sizeof(w32));   cntR++; cntRpad++; cntRrow++;   cntRdig++;//Printf(" WR: %x isDig: %d",w32,isDigThere[iddl][row][dil][pad]);
351               }
352             }//ddl 
353           }//isDig
354         }//pad
355         WriteEoE(ddlL,row,dil,cntLpad); cntL++;  cntLrow++;    cntLeoe++;                                 //molnarl: write EoE markers
356         WriteEoE(ddlR,row,dil,cntRpad); cntR++;  cntRrow++;    cntReoe++;
357       }//dil
358       if(row%8==0){                                               
359         WriteSegMarker(ddlL,row); cntL++;  cntLseg++;
360         WriteSegMarker(ddlR,row); cntR++;  cntRseg++;
361       }
362       posL=ddlL->Tellp();   ddlL->Seekp(posLmarker);    WriteRowMarker(ddlL,(UInt_t)(cntLrow-1)); ddlL->Seekp(posL);      //find the marker position write and  go back to the actual position to continue writing                    
363       posR=ddlR->Tellp();   ddlR->Seekp(posRmarker);    WriteRowMarker(ddlR,(UInt_t)(cntRrow-1)); ddlR->Seekp(posR);                           
364     }//row
365     header.fSize=sizeof(header)+cntL*sizeof(w32); ddlL->Seekp(0); ddlL->WriteBuffer((char*)&header,sizeof(header)); delete ddlL; //rewrite header with size set to
366     header.fSize=sizeof(header)+cntR*sizeof(w32); ddlR->Seekp(0); ddlR->WriteBuffer((char*)&header,sizeof(header)); delete ddlR; //number of bytes and close file
367     
368     //Printf("In Ch %d # digits written to LDD %d RDDL %d",iCh,cntLdig,cntRdig);
369     
370   }//chambers loop
371 }//WriteRaw()
372 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
373
374 #endif