]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/readSSDOCDBEntry.C
New task to check TPC-ITS track prolongation eff with cosmics
[u/mrichter/AliRoot.git] / ITS / readSSDOCDBEntry.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include "TString.h"
3 #include "TH1.h"
4 #include "TH2.h"
5 #include "TFile.h"
6 #include "TObjArray.h"
7 #include "TObject.h"
8 #include "TCanvas.h"
9 #include "TStyle.h"
10 #include "Riostream.h"
11
12 #include "AliITSgeomTGeo.h"
13 #include "AliITSGainSSDv2.h"
14 #include "AliITSBadChannelsSSDv2.h"
15 #include "AliITSNoiseSSDv2.h"
16 #include "AliITSGainSSD.h"
17 #include "AliITSBadChannelsSSD.h"
18 #include "AliITSNoiseSSD.h"
19 #include "AliCDBManager.h"
20 #include "AliCDBEntry.h"
21 #endif
22
23 /*  $Id$    */
24
25
26
27 //====================================================================//
28 void Noise(AliCDBManager * man, Int_t runNumber);
29 void BadChannelMap(AliCDBManager * man);
30 void GainCalibration(AliCDBManager * man);
31 void ReadOldSSDNoise(TObjArray *array, AliITSNoiseSSDv2 *noiseSSD);
32 void ReadOldSSDBadChannels(TObjArray *array, AliITSBadChannelsSSDv2 *badChannelsSSD);
33 void ReadOldSSDGain(TObjArray *array, AliITSGainSSDv2 *gainSSD);
34 //====================================================================//
35
36 //_____________________________________________________________________//
37 void readSSDOCDBEntry(const char* type = "alien", Int_t runNumber = 0) {
38   //This macro allows to visualize the bad channels in the OCDB
39   //The type can be either "local" or "alien" (where the OCDB file comes from)
40   //The run nmber is the pedestal one
41   gStyle->SetPalette(1,0);
42   TString gType = type;
43   
44   AliCDBManager * man = AliCDBManager::Instance();
45   
46   if(gType == "alien") {
47     man->SetDefaultStorage("alien://folder=/alice/data/2009/OCDB/");
48   }
49   else if(gType == "local") 
50     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
51   else {
52     cout<<"Allowed types: local or alien!!!"<<endl;
53     abort();
54   }
55   
56   man->SetRun(runNumber);
57     
58   Noise(man,runNumber);
59   BadChannelMap(man);
60   GainCalibration(man);
61 }
62
63 //_____________________________________________________________________//
64 void Noise(AliCDBManager * man, Int_t runNumber) {
65   const Int_t fgkSSDMODULES = 1698;
66   const Int_t fgkSSDSTRIPSPERMODULE = 1536;
67   static const Int_t fgkDefaultNStripsSSD = 768;
68
69   //noise histograms
70   Int_t fLayer = 0,fLadder = 0, fModule = 0;
71   Int_t fHistCounter = 0;
72   TString fTitle;
73   TObjArray *array = new TObjArray();
74   TH1D *hNoiseModule[fgkSSDMODULES];   
75   for(Int_t i = 500; i < fgkSSDMODULES + 500; i++) {
76     AliITSgeomTGeo::GetModuleId(i,fLayer,fLadder,fModule);
77     fTitle = "SSD_Noise_Layer"; fTitle += fLayer;
78     fTitle += "_Ladder"; fTitle += fLadder;
79     fTitle += "_Module"; fTitle += fModule;
80     
81     hNoiseModule[fHistCounter] = new TH1D(fTitle.Data(),fTitle.Data(),1540,0,1540);
82     hNoiseModule[fHistCounter]->GetXaxis()->SetTitleColor(1);
83     hNoiseModule[fHistCounter]->GetXaxis()->SetTitle("Strip number");
84     hNoiseModule[fHistCounter]->GetYaxis()->SetTitle("Noise");
85     array->AddLast(hNoiseModule[fHistCounter]);
86     fHistCounter += 1;
87   }
88   
89   AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();
90   AliCDBEntry *entryNoiseSSD = man->Get("ITS/Calib/NoiseSSD");
91   TObject *empty = (TObject *)entryNoiseSSD->GetObject();
92   TString objectname = empty->GetName();
93   if(objectname=="TObjArray") {
94     TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
95     ReadOldSSDNoise(noiseSSDOld, noiseSSD);
96   }
97   else if(objectname=="AliITSNoiseSSDv2") {
98     cout<<"Reading the new format of the calibration file..."<<endl;
99     noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
100   }
101   
102   Double_t noise = 0.0;
103   for (Int_t i = 0; i < fgkSSDMODULES; i++) {
104     //cout<<"Noise for module: "<<i+1<<endl;
105     for(Int_t j = 0; j < fgkDefaultNStripsSSD; j++) {
106       noise = noiseSSD->GetNoiseP(i,j);
107       hNoiseModule[i]->SetBinContent(j+1,noise);
108       noise = noiseSSD->GetNoiseN(i,j);
109       hNoiseModule[i]->SetBinContent(fgkSSDSTRIPSPERMODULE-j,noise);
110     }//loop over strips
111   }//loop over modules
112   
113   TString output = "noiseSSD."; output += runNumber; output += ".root";
114   TFile *f = TFile::Open(output.Data(),"recreate");
115   array->Write();
116   f->Close();
117 }
118
119 //_____________________________________________________________________//
120 void BadChannelMap(AliCDBManager * man) {
121   const Int_t fgkSSDMODULES = 1698;
122   static const Int_t fgkDefaultNStripsSSD = 768;
123
124   //_____________________________________________________________________________//
125   TH2F *fHistPSideBadChannelMapLayer5 = new TH2F("fHistPSideBadChannelMapLayer5",
126                                                  "Layer 5;N_{module};N_{ladder}",
127                                                  22,1,23,
128                                                  34,500,534);
129   fHistPSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
130   fHistPSideBadChannelMapLayer5->GetZaxis()->SetRangeUser(0.,100.);
131   fHistPSideBadChannelMapLayer5->SetStats(kFALSE);
132   fHistPSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
133   fHistPSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
134   fHistPSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
135   fHistPSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
136   fHistPSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
137   fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
138   fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
139   TH2F *fHistNSideBadChannelMapLayer5 = new TH2F("fHistNSideBadChannelMapLayer5",
140                                                  "Layer 5;N_{module};N_{ladder}",
141                                                  22,1,23,
142                                                  34,500,534);
143   fHistNSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
144   fHistNSideBadChannelMapLayer5->GetZaxis()->SetRangeUser(0.,100.);
145   fHistNSideBadChannelMapLayer5->SetStats(kFALSE);
146   fHistNSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
147   fHistNSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
148   fHistNSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
149   fHistNSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
150   fHistNSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
151   fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
152   fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
153
154   TH2F *fHistPSideBadChannelMapLayer6 = new TH2F("fHistPSideBadChannelMapLayer6",
155                                             "Layer 6;N_{module};N_{ladder}",
156                                             25,1,26,
157                                             38,600,638);
158   fHistPSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
159   fHistPSideBadChannelMapLayer6->GetZaxis()->SetRangeUser(0.,100.);
160   fHistPSideBadChannelMapLayer6->SetStats(kFALSE);
161   fHistPSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
162   fHistPSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
163   fHistPSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
164   fHistPSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
165   fHistPSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
166   fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
167   fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
168   TH2F *fHistNSideBadChannelMapLayer6 = new TH2F("fHistNSideBadChannelMapLayer6",
169                                             "Layer 6;N_{module};N_{ladder}",
170                                             25,1,26,
171                                             38,600,638);
172   fHistNSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
173   fHistNSideBadChannelMapLayer6->GetZaxis()->SetRangeUser(0.,100.);
174   fHistNSideBadChannelMapLayer6->SetStats(kFALSE);
175   fHistNSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
176   fHistNSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
177   fHistNSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
178   fHistNSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
179   fHistNSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
180   fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
181   fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
182   //_____________________________________________________________________________//
183   
184   //_____________________________________________________________________________//
185   AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
186   AliCDBEntry *entryBadChannelsSSD = man->Get("ITS/Calib/BadChannelsSSD");
187   TObject *empty = (TObject *)entryBadChannelsSSD->GetObject();
188   TString objectname = empty->GetName();
189   if(objectname=="TObjArray") {
190     TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
191     ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
192   }
193   else if(objectname=="AliITSBadChannelsSSDv2") {
194     cout<<"Reading the new format of the calibration file..."<<endl;
195     badChannelsSSD = (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
196   }
197   //_____________________________________________________________________________//
198
199   //_____________________________________________________________________________//
200   Int_t nPSideChannelsTotal = 0, nNSideChannelsTotal = 0;
201   Int_t nBadPSideChannelsTotal = 0, nBadNSideChannelsTotal = 0;
202   Int_t nBadPSideChannels = 0, nBadNSideChannels = 0;
203   Int_t layer = 0, ladder = 0, module = 0;
204   Int_t nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
205   Int_t nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
206   //_____________________________________________________________________________//
207
208   for(Int_t i = 0; i < fgkSSDMODULES; i++) {
209     //for(Int_t i = 0; i < 1; i++) {
210     AliITSgeomTGeo::GetModuleId(i+500,layer,ladder,module);
211     nBadPSideChannels = 0, nBadNSideChannels = 0;
212     nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
213     nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
214
215     Int_t badChannel = 0;
216     for(Int_t j = 0; j < fgkDefaultNStripsSSD; j++) {
217       badChannel = (Int_t)(badChannelsSSD->GetBadChannelP(i,j));
218       //cout<<"Module: "<<i+500<< " Strip: "<<j<<" - "<<badChannel<<endl;
219       if(badChannel != 0) {
220         if(layer == 5)
221           nPSideChannelsLayer5 += 1;
222         if(layer == 6)
223           nPSideChannelsLayer6 += 1;
224         nBadPSideChannels += 1;
225       }
226       badChannel = (Int_t)(badChannelsSSD->GetBadChannelN(i,j));
227       //cout<<"Module: "<<i+500<< " Strip: "<<fgkDefaultNStripsSSD+j+1<<" - "<<badChannel<<endl;
228       if(badChannel != 0) {
229         if(layer == 5)                                                    
230           nNSideChannelsLayer5 += 1;
231         if(layer == 6)
232           nNSideChannelsLayer6 += 1;
233         nBadNSideChannels += 1;
234       }
235     }
236     if(layer == 5) {
237       if(nPSideChannelsLayer5 > 0)
238         fHistPSideBadChannelMapLayer5->Fill(module,499+ladder,
239                                             100.*nPSideChannelsLayer5/fgkDefaultNStripsSSD);
240       else fHistPSideBadChannelMapLayer5->Fill(module,499+ladder,0.0001);
241       if(nNSideChannelsLayer5 > 0)
242         fHistNSideBadChannelMapLayer5->Fill(module,499+ladder,
243                                             100.*nNSideChannelsLayer5/fgkDefaultNStripsSSD);
244       else fHistNSideBadChannelMapLayer5->Fill(module,499+ladder,0.0001);
245     }//layer 5
246     if(layer == 6) {
247       if(nPSideChannelsLayer6 > 0) 
248         fHistPSideBadChannelMapLayer6->Fill(module,599+ladder,
249                                             100.*nPSideChannelsLayer6/fgkDefaultNStripsSSD);
250       else fHistPSideBadChannelMapLayer6->Fill(module,599+ladder,0.0001);
251       if(nNSideChannelsLayer6 > 0) 
252         fHistNSideBadChannelMapLayer6->Fill(module,599+ladder,
253                                             100.*nNSideChannelsLayer6/fgkDefaultNStripsSSD);
254       else fHistNSideBadChannelMapLayer6->Fill(module,599+ladder,0.0001);
255     }//layer 6
256       
257     nBadPSideChannelsTotal += nBadPSideChannels;
258     nBadNSideChannelsTotal += nBadNSideChannels;
259     nPSideChannelsTotal += fgkDefaultNStripsSSD;
260     nNSideChannelsTotal += fgkDefaultNStripsSSD;
261   }
262
263   cout<<"================================="<<endl;
264   cout<<"Bad p-Side channels: "<<100.*nBadPSideChannelsTotal/nPSideChannelsTotal<<endl;
265   cout<<"Bad n-Side channels: "<<100.*nBadNSideChannelsTotal/nNSideChannelsTotal<<endl;
266   cout<<"================================="<<endl;
267
268   TCanvas *cBadChannel = new TCanvas("cBadChannel",
269                                      "Bad channel list",0,0,900,900);
270   cBadChannel->SetHighLightColor(10); cBadChannel->SetFillColor(10); 
271   cBadChannel->Divide(2,2);
272
273   cBadChannel->cd(1)->SetBottomMargin(.2); 
274   cBadChannel->cd(1)->SetLeftMargin(.15);
275   cBadChannel->cd(1)->SetRightMargin(.2);
276   cBadChannel->cd(1)->SetGridx(); cBadChannel->cd(1)->SetGridy();
277   cBadChannel->cd(1); fHistPSideBadChannelMapLayer5->Draw("colz"); 
278   cBadChannel->cd(2)->SetBottomMargin(.2); 
279   cBadChannel->cd(2)->SetLeftMargin(.15);
280   cBadChannel->cd(2)->SetRightMargin(.2);
281   cBadChannel->cd(2)->SetGridx(); cBadChannel->cd(2)->SetGridy();
282   cBadChannel->cd(2); fHistPSideBadChannelMapLayer6->Draw("colz");
283   cBadChannel->cd(3)->SetBottomMargin(.2); 
284   cBadChannel->cd(3)->SetLeftMargin(.15);
285   cBadChannel->cd(3)->SetRightMargin(.2);
286   cBadChannel->cd(3)->SetGridx(); cBadChannel->cd(3)->SetGridy();
287   cBadChannel->cd(3); fHistNSideBadChannelMapLayer5->Draw("colz"); 
288   cBadChannel->cd(4)->SetBottomMargin(.2); 
289   cBadChannel->cd(4)->SetLeftMargin(.15);
290   cBadChannel->cd(4)->SetRightMargin(.2);
291   cBadChannel->cd(4)->SetGridx(); cBadChannel->cd(4)->SetGridy();
292   cBadChannel->cd(4); fHistNSideBadChannelMapLayer6->Draw("colz");
293   cBadChannel->SaveAs("Run-BadChannels.gif");
294 }
295
296 //_____________________________________________________________________//
297 void GainCalibration(AliCDBManager * man) {
298   const Int_t fgkSSDMODULES = 1698;
299   static const Int_t fgkDefaultNStripsSSD = 768;
300
301   TH2F *fHistGainMapLayer5 = new TH2F("fHistGainMapLayer5",
302                                       "Layer 5;N_{strip};N_{module}",
303                                       1537,0,1537,
304                                       750,499,1249);
305   fHistGainMapLayer5->GetXaxis()->SetTitleColor(1);
306   fHistGainMapLayer5->SetStats(kFALSE);
307   fHistGainMapLayer5->GetYaxis()->SetTitleOffset(1.8);
308   TH2F *fHistGainMapLayer6 = new TH2F("fHistGainMapLayer6",
309                                       "Layer 6;N_{strip};N_{module}",
310                                       1537,0,1537,
311                                       952,1249,2199);
312   fHistGainMapLayer6->GetXaxis()->SetTitleColor(1);
313   fHistGainMapLayer6->SetStats(kFALSE);
314   fHistGainMapLayer6->GetYaxis()->SetTitleOffset(1.8);
315
316   AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
317   AliCDBEntry *entryGainSSD = man->Get("ITS/Calib/GainSSD");
318   TObject *empty = (TObject *)entryGainSSD->GetObject();
319   TString objectname = empty->GetName();
320   if(objectname=="Gain") {
321     TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
322     ReadOldSSDGain(gainSSDOld, gainSSD);
323   }
324   else if(objectname=="AliITSGainSSDv2") {
325     cout<<"Reading the new format of the calibration file..."<<endl;
326     gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
327   }
328   
329   Int_t layer = 0, ladder = 0, module = 0;
330   Double_t gain = 0;
331   for (Int_t i = 0; i < fgkSSDMODULES; i++) {
332   //for (Int_t i = 0; i < 1; i++) {
333     AliITSgeomTGeo::GetModuleId(i+500,layer,ladder,module);
334     for(Int_t j = 0; j < fgkDefaultNStripsSSD; j++) {
335       gain = gainSSD->GetGainP(i,j);
336       //cout<<"GainP: "<<gain<<endl;
337       if(layer == 5)
338         fHistGainMapLayer5->Fill(j,i+500,gain);
339       if(layer == 6)
340         fHistGainMapLayer6->Fill(j,i+500,gain);
341       
342       gain = gainSSD->GetGainN(i,j);
343       //cout<<"GainN: "<<gain<<endl;
344       if(layer == 5)
345         fHistGainMapLayer5->Fill(fgkDefaultNStripsSSD+j,i+500,gain);
346       if(layer == 6)
347         fHistGainMapLayer6->Fill(fgkDefaultNStripsSSD+j,i+500,gain);
348     }//strip loop
349   }//module loop
350
351   TCanvas *cGain = new TCanvas("cGain","Gain calibration map",0,300,600,300);
352   cGain->SetHighLightColor(10); cGain->SetFillColor(10); cGain->Divide(2,1);
353   
354   cGain->cd(1)->SetBottomMargin(.2); cGain->cd(1)->SetLeftMargin(.15);
355   cGain->cd(1); fHistGainMapLayer5->Draw("colz");
356   cGain->cd(2)->SetBottomMargin(.2); cGain->cd(2)->SetLeftMargin(.15);
357   cGain->cd(2); fHistGainMapLayer6->Draw("colz");
358 }
359
360 //_____________________________________________________________________//
361 void ReadOldSSDNoise(TObjArray *array, 
362                      AliITSNoiseSSDv2 *noiseSSD) {
363   const Int_t fgkSSDSTRIPSPERMODULE = 1536;
364   const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
365
366   Int_t fNMod = array->GetEntries();
367   cout<<"Converting old calibration object for noise..."<<endl;
368
369   //NOISE
370   Double_t noise = 0.0;
371   for (Int_t iModule = 0; iModule < fNMod; iModule++) {
372     AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
373     for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
374       noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
375       if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
376         noiseSSD->AddNoiseP(iModule,iStrip,noise);
377       if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
378         noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
379     }//loop over strips
380   }//loop over modules      
381 }
382
383 //_____________________________________________________________________//
384 void ReadOldSSDBadChannels(TObjArray *array, 
385                            AliITSBadChannelsSSDv2 *badChannelsSSD) {
386   Int_t fNMod = array->GetEntries();
387   cout<<"Converting old calibration object for bad channels..."<<endl;
388   for (Int_t iModule = 0; iModule < fNMod; iModule++) {
389     //for (Int_t iModule = 0; iModule < 1; iModule++) {
390     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
391     TArrayI arrayPSide = bad->GetBadPChannelsList();
392     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++) 
393       badChannelsSSD->AddBadChannelP(iModule,
394                                      iPCounter,
395                                      (Char_t)arrayPSide.At(iPCounter));
396         
397     TArrayI arrayNSide = bad->GetBadNChannelsList();
398     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++) 
399       badChannelsSSD->AddBadChannelN(iModule,
400                                      iNCounter,
401                                      (Char_t)arrayNSide.At(iNCounter));
402     
403   }//loop over modules      
404 }
405
406 //_____________________________________________________________________//
407 void ReadOldSSDGain(TObjArray *array, 
408                     AliITSGainSSDv2 *gainSSD) {
409   Int_t fNMod = array->GetEntries();
410   cout<<"Converting old calibration object for gain..."<<endl;
411
412   //GAIN
413   for (Int_t iModule = 0; iModule < fNMod; iModule++) {
414     AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
415     TArrayF arrayPSide = gainModule->GetGainP();
416     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
417       gainSSD->AddGainP(iModule,
418                         iPCounter,
419                         arrayPSide.At(iPCounter));
420     TArrayF arrayNSide = gainModule->GetGainN();
421     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
422       gainSSD->AddGainN(iModule,
423                         iNCounter,
424                         arrayNSide.At(iNCounter));
425   }//loop over modules 
426 }