]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliCaloCalibPedestal.cxx
remove unused AliEMCALGeometryOfflineTrd1
[u/mrichter/AliRoot.git] / EMCAL / AliCaloCalibPedestal.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 /* $Id$ */
16
17 //________________________________________________________________________
18 //
19 // A help class for monitoring and calibration tools: MOOD, AMORE etc.,
20 // It can be created and used a la (ctor):
21 /*
22   //Create the object for making the histograms
23   fPedestals = new AliCaloCalibPedestal( fDetType );
24   // AliCaloCalibPedestal knows how many modules we have for PHOS or EMCAL
25   fNumModules = fPedestals->GetModules();
26 */
27 // fed an event:
28 //  fPedestals->ProcessEvent(fCaloRawStream);
29 // asked to draw histograms:
30 //  fPedestals->GetDeadMap(i)->Draw("col");
31 // or
32 //  fPedestals->GetPeakProfileHighGainRatio((i < fNumModules) ? i : fVisibleModule)->Draw("colz");
33 // etc.
34 // The pseudo-code examples above were from the first implementation in MOOD (summer 2007).
35 //________________________________________________________________________
36
37 #include "TH1.h"
38 #include "TFile.h"
39
40 #include <fstream>
41 #include <iostream>
42
43 #include "AliCaloRawStream.h"
44
45 //The include file
46 #include "AliCaloCalibPedestal.h"
47
48 ClassImp(AliCaloCalibPedestal)
49
50 using namespace std;
51
52 // ctor; initialize everything in order to avoid compiler warnings
53 AliCaloCalibPedestal::AliCaloCalibPedestal(kDetType detectorType) :  
54   TObject(),
55   fPedestalLowGain(),
56   fPedestalHighGain(),
57   fPeakMinusPedLowGain(),
58   fPeakMinusPedHighGain(),
59   fPedestalLowGainDiff(),
60   fPedestalHighGainDiff(),
61   fPeakMinusPedLowGainDiff(),
62   fPeakMinusPedHighGainDiff(),
63   fPedestalLowGainRatio(),
64   fPedestalHighGainRatio(),
65   fPeakMinusPedLowGainRatio(),
66   fPeakMinusPedHighGainRatio(),
67   fDeadMap(),
68   fDeadTowers(0),
69   fNewDeadTowers(0),
70   fResurrectedTowers(0),
71   fReference(0),
72   fDetType(kNone),
73   fColumns(0),
74   fRows(0),
75   fModules(0),
76   fRunNumber(-1)
77 {
78   //Default constructor. First we set the detector-type related constants.
79   if (detectorType == kPhos) {
80     fColumns = fgkPhosCols;
81     fRows = fgkPhosRows;
82     fModules = fgkPhosModules;
83   } 
84   else {
85     //We'll just trust the enum to keep everything in line, so that if detectorType
86     //isn't kPhos then it is kEmCal. Note, however, that this is not necessarily the
87     //case, if someone intentionally gives another number
88     fColumns = fgkEmCalCols;
89     fRows = fgkEmCalRows;
90     fModules = fgkEmCalModules;
91   } 
92   fDetType = detectorType;
93  
94   //Then, loop for the requested number of modules
95   TString title, name;
96   for (int i = 0; i < fModules; i++) {
97     //Pedestals, low gain
98     name = "hPedlowgain";
99     name += i;
100     title = "Pedestals, low gain, module ";
101     title += i; 
102     fPedestalLowGain.Add(new TProfile2D(name, title,
103                                         fColumns, 0.0, fColumns, 
104                                         fRows, -fRows, 0.0));
105   
106     //Pedestals, high gain
107     name = "hPedhighgain";
108     name += i;
109     title = "Pedestals, high gain, module ";
110     title += i; 
111     fPedestalHighGain.Add(new TProfile2D(name, title,
112                                          fColumns, 0.0, fColumns, 
113                                          fRows, -fRows, 0.0));
114   
115     //Peak-Pedestals, low gain
116     name = "hPeakMinusPedlowgain";
117     name += i;
118     title = "Peak-Pedestal, low gain, module ";
119     title += i; 
120     fPeakMinusPedLowGain.Add(new TProfile2D(name, title,
121                                             fColumns, 0.0, fColumns, 
122                                             fRows, -fRows, 0.0));
123   
124     //Peak-Pedestals, high gain
125     name = "hPeakMinusPedhighgain";
126     name += i;
127     title = "Peak-Pedestal, high gain, module ";
128     title += i; 
129     fPeakMinusPedHighGain.Add(new TProfile2D(name, title,
130                                              fColumns, 0.0, fColumns, 
131                                              fRows, -fRows, 0.0));
132   
133     name = "hDeadMap";
134     name += i;
135     title = "Dead map, module ";
136     title += i;
137     fDeadMap.Add(new TH2D(name, title, fColumns, 0.0, fColumns, 
138                           fRows, -fRows, 0.0));
139   
140   }//end for nModules create the histograms
141  
142   //Compress the arrays, in order to remove the empty objects (a 16 slot array is created by default)
143   fPedestalLowGain.Compress();
144   fPedestalHighGain.Compress();
145   fPeakMinusPedLowGain.Compress();
146   fPeakMinusPedHighGain.Compress();
147   fDeadMap.Compress();
148   //Make them the owners of the profiles, so we don't need to care about deleting them
149   //fPedestalLowGain.SetOwner();
150   //fPedestalHighGain.SetOwner();
151   //fPeakMinusPedLowGain.SetOwner();
152   //fPeakMinusPedHighGain.SetOwner();
153   
154 }
155
156 // dtor
157 //_____________________________________________________________________
158 AliCaloCalibPedestal::~AliCaloCalibPedestal()
159 {
160   if (fReference) delete fReference;//Delete the reference object, if it has been loaded
161   //TObjArray will delete the histos/profiles when it is deleted.
162 }
163
164 // copy ctor
165 //_____________________________________________________________________
166 AliCaloCalibPedestal::AliCaloCalibPedestal(const AliCaloCalibPedestal &ped) :
167   TObject(ped),
168   fPedestalLowGain(),
169   fPedestalHighGain(),
170   fPeakMinusPedLowGain(),
171   fPeakMinusPedHighGain(),
172   fPedestalLowGainDiff(),
173   fPedestalHighGainDiff(),
174   fPeakMinusPedLowGainDiff(),
175   fPeakMinusPedHighGainDiff(),
176   fPedestalLowGainRatio(),
177   fPedestalHighGainRatio(),
178   fPeakMinusPedLowGainRatio(),
179   fPeakMinusPedHighGainRatio(),
180   fDeadMap(),
181   fDeadTowers(ped.GetDeadTowerCount()),
182   fNewDeadTowers(ped.GetDeadTowerNew()),
183   fResurrectedTowers(ped.GetDeadTowerResurrected()),
184   fReference( 0 ), //! note that we do not try to copy the reference info here
185   fDetType(ped.GetDetectorType()),
186   fColumns(ped.GetColumns()),
187   fRows(ped.GetRows()),
188   fModules(ped.GetModules()),
189   fRunNumber(ped.GetRunNumber())
190 {
191   // Then the ObjArray ones; we add the histograms rather than trying TObjArray = assignment
192   //DS: this has not really been tested yet..
193   for (int i = 0; i < fModules; i++) {
194     fPedestalLowGain.Add( ped.GetPedProfileLowGain(i) );
195     fPedestalHighGain.Add( ped.GetPedProfileHighGain(i) );
196     fPeakMinusPedLowGain.Add( ped.GetPeakProfileLowGain(i) );
197     fPeakMinusPedHighGain.Add( ped.GetPeakProfileHighGain(i) );
198
199     fDeadMap.Add( ped.GetDeadMap(i) );  
200   }//end for nModules 
201  
202   //Compress the arrays, in order to remove the empty objects (a 16 slot array is created by default)
203   fPedestalLowGain.Compress();
204   fPedestalHighGain.Compress();
205   fPeakMinusPedLowGain.Compress();
206   fPeakMinusPedHighGain.Compress();
207   fDeadMap.Compress();
208 }
209
210 // assignment operator; use copy ctor to make life easy..
211 //_____________________________________________________________________
212 AliCaloCalibPedestal& AliCaloCalibPedestal::operator = (const AliCaloCalibPedestal &source)
213 {
214   // assignment operator; use copy ctor
215   if (&source == this) return *this;
216
217   new (this) AliCaloCalibPedestal(source);
218   return *this;
219 }
220
221 //_____________________________________________________________________
222 void AliCaloCalibPedestal::Reset()
223 {
224   // Reset all arrays/histograms
225   for (int i = 0; i < fModules; i++) {
226     GetPedProfileLowGain(i)->Reset();
227     GetPedProfileHighGain(i)->Reset();
228     GetPeakProfileLowGain(i)->Reset();
229     GetPeakProfileHighGain(i)->Reset();
230     GetDeadMap(i)->Reset();
231     
232     if (!fPedestalLowGainDiff.IsEmpty()) {
233       //This means that the comparison profiles have been created.
234   
235       GetPedProfileLowGainDiff(i)->Reset();
236       GetPedProfileHighGainDiff(i)->Reset();
237       GetPeakProfileLowGainDiff(i)->Reset();
238       GetPeakProfileHighGainDiff(i)->Reset();
239       
240       GetPedProfileLowGainRatio(i)->Reset();
241       GetPedProfileHighGainRatio(i)->Reset();
242       GetPeakProfileLowGainRatio(i)->Reset();
243       GetPeakProfileHighGainRatio(i)->Reset();
244     }
245   }
246   fDeadTowers = 0;
247   fNewDeadTowers = 0;
248   fResurrectedTowers = 0;
249  
250   //To think about: should fReference be deleted too?... let's not do it this time, at least...
251 }
252
253 //_____________________________________________________________________
254 Bool_t AliCaloCalibPedestal::ProcessEvent(AliCaloRawStream *in)
255
256   // Method to process=analyze one event in the data stream
257   if (!in) return kFALSE; //Return right away if there's a null pointer
258   
259   in->SetOldRCUFormat(kTRUE);
260   
261   int sample, i = 0; //The sample temp, and the sample number in current event.
262   int max = fgkSampleMin, min = fgkSampleMax;//Use these for picking the pedestal
263   int gain = 0;
264   
265   while (in->Next()) {
266     sample = in->GetSignal(); //Get the adc signal
267     if (sample < min) min = sample;
268     if (sample > max) max = sample;
269     i++;
270     if ( i >= in->GetTimeLength()) {
271       //If we're here then we're done with this tower
272       gain = 1 - in->IsLowGain();
273       
274       int arrayPos = in->GetModule(); //The modules are numbered starting from 0
275       if (arrayPos >= fModules) {
276         //TODO: return an error message, if appopriate (perhaps if debug>0?)
277         return kFALSE;
278       } 
279     
280       //Debug
281       if (arrayPos < 0 || arrayPos >= fModules) {
282         printf("Oh no: arrayPos = %i.\n", arrayPos); 
283       }
284
285       //NOTE: coordinates are (column, row) for the profiles
286       if (gain == 0) {
287         //fill the low gain histograms
288         ((TProfile2D*)fPedestalLowGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, min);
289         ((TProfile2D*)fPeakMinusPedLowGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, max - min);
290       } 
291       else {//fill the high gain ones
292         ((TProfile2D*)fPedestalHighGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, min);
293         ((TProfile2D*)fPeakMinusPedHighGain[arrayPos])->Fill(in->GetColumn(), -in->GetRow() - 1, max - min);
294       }//end if gain
295       
296       max = fgkSampleMin; min = fgkSampleMax;
297       i = 0;
298     
299     }//End if end of tower
300    
301   }//end while, of stream
302   
303   return kTRUE;
304 }
305
306 //_____________________________________________________________________
307 Bool_t AliCaloCalibPedestal::SaveHistograms(TString fileName, Bool_t saveEmptyHistos)
308 {
309   //Saves all the histograms (or profiles, to be accurate) to the designated file
310   
311   TFile destFile(fileName, "recreate");
312   
313   if (destFile.IsZombie()) {
314     return kFALSE;
315   }
316   
317   destFile.cd();
318   
319   for (int i = 0; i < fModules; i++) {
320     if( ((TProfile2D *)fPeakMinusPedLowGain[i])->GetEntries() || saveEmptyHistos) {
321       fPeakMinusPedLowGain[i]->Write();
322     }
323     if( ((TProfile2D *)fPeakMinusPedHighGain[i])->GetEntries() || saveEmptyHistos) { 
324       fPeakMinusPedHighGain[i]->Write();
325     }
326     if( ((TProfile2D *)fPedestalLowGain[i])->GetEntries() || saveEmptyHistos) {
327       fPedestalLowGain[i]->Write();
328     }
329     if( ((TProfile2D *)fPedestalHighGain[i])->GetEntries() || saveEmptyHistos) {
330       fPedestalHighGain[i]->Write();
331     }
332   } 
333   
334   destFile.Close();
335   
336   return kTRUE;
337 }
338
339 //_____________________________________________________________________
340 Bool_t AliCaloCalibPedestal::LoadReferenceCalib(TString fileName, TString objectName)
341 {
342   
343   //Make sure that the histograms created when loading the object are not destroyed as the file object is destroyed
344   TH1::AddDirectory(kFALSE);
345   
346   TFile *sourceFile = new TFile(fileName);
347   if (sourceFile->IsZombie()) {
348     return kFALSE;//We couldn't load the reference
349   }
350
351   if (fReference) delete fReference;//Delete the reference object, if it already exists
352   fReference = 0;
353   
354   fReference = (AliCaloCalibPedestal*)sourceFile->Get(objectName);
355  
356   if (!fReference || !(fReference->InheritsFrom(AliCaloCalibPedestal::Class())) || (fReference->GetDetectorType() != fDetType)) {
357     if (fReference) delete fReference;//Delete the object, in case we had an object of the wrong type
358     fReference = 0;
359     return kFALSE;
360   }
361         
362   delete sourceFile;
363  
364   //Reset the histogram ownership behaviour. NOTE: a better workaround would be good, since this may accidentally set AddDirectory to true, even
365   //if we are called by someone who has set it to false...
366   TH1::AddDirectory(kTRUE);
367  
368   return kTRUE;//We succesfully loaded the object
369 }
370
371 //_____________________________________________________________________
372 void AliCaloCalibPedestal::ValidateComparisonProfiles()
373 {
374   //Make sure the comparison histos exist
375   if (!fPedestalLowGainDiff.IsEmpty()) return; //The profiles already exist. We just check one, because they're all created at
376   //the same time
377                                                 
378                                                 
379   //Then, loop for the requested number of modules
380   TString title, name;
381   for (int i = 0; i < fModules; i++) {
382     //Pedestals, low gain
383     name = "hPedlowgainDiff";
384     name += i;
385     title = "Pedestals difference, low gain, module ";
386     title += i; 
387     fPedestalLowGainDiff.Add(new TProfile2D(name, title,
388                                             fColumns, 0.0, fColumns, 
389                                             fRows, -fRows, 0.0));
390   
391     //Pedestals, high gain
392     name = "hPedhighgainDiff";
393     name += i;
394     title = "Pedestals difference, high gain, module ";
395     title += i; 
396     fPedestalHighGainDiff.Add(new TProfile2D(name, title,
397                                              fColumns, 0.0, fColumns, 
398                                              fRows, -fRows, 0.0));
399   
400     //Peak-Pedestals, low gain
401     name = "hPeakMinusPedlowgainDiff";
402     name += i;
403     title = "Peak-Pedestal difference, low gain, module ";
404     title += i; 
405     fPeakMinusPedLowGainDiff.Add(new TProfile2D(name, title,
406                                                 fColumns, 0.0, fColumns, 
407                                                 fRows, -fRows, 0.0));
408   
409     //Peak-Pedestals, high gain
410     name = "hPeakMinusPedhighgainDiff";
411     name += i;
412     title = "Peak-Pedestal difference, high gain, module ";
413     title += i; 
414     fPeakMinusPedHighGainDiff.Add(new TProfile2D(name, title,
415                                                  fColumns, 0.0, fColumns, 
416                                                  fRows, -fRows, 0.0));
417   
418     //Pedestals, low gain
419     name = "hPedlowgainRatio";
420     name += i;
421     title = "Pedestals ratio, low gain, module ";
422     title += i; 
423     fPedestalLowGainRatio.Add(new TProfile2D(name, title,
424                                              fColumns, 0.0, fColumns, 
425                                              fRows, -fRows, 0.0));
426   
427     //Pedestals, high gain
428     name = "hPedhighgainRatio";
429     name += i;
430     title = "Pedestals ratio, high gain, module ";
431     title += i; 
432     fPedestalHighGainRatio.Add(new TProfile2D(name, title,
433                                               fColumns, 0.0, fColumns, 
434                                               fRows, -fRows, 0.0));
435   
436     //Peak-Pedestals, low gain
437     name = "hPeakMinusPedlowgainRatio";
438     name += i;
439     title = "Peak-Pedestal ratio, low gain, module ";
440     title += i; 
441     fPeakMinusPedLowGainRatio.Add(new TProfile2D(name, title,
442                                                  fColumns, 0.0, fColumns, 
443                                                  fRows, -fRows, 0.0));
444   
445     //Peak-Pedestals, high gain
446     name = "hPeakMinusPedhighgainRatio";
447     name += i;
448     title = "Peak-Pedestal ratio, high gain, module ";
449     title += i; 
450     fPeakMinusPedHighGainRatio.Add(new TProfile2D(name, title,
451                                                   fColumns, 0.0, fColumns, 
452                                                   fRows, -fRows, 0.0));
453     
454   }//end for nModules create the histograms
455 }
456
457 //_____________________________________________________________________
458 void AliCaloCalibPedestal::ComputeDiffAndRatio()
459 {
460   // calculate differences and ratios relative to a reference
461   ValidateComparisonProfiles();//Make sure the comparison histos exist
462  
463   if (!fReference) {
464     return;//Return if the reference object isn't loaded
465   }
466
467   for (int i = 0; i < fModules; i++) {
468     //Compute the ratio of the histograms
469     
470     ((TProfile2D*)fPedestalLowGainRatio[i])->Divide(GetPedProfileLowGain(i), fReference->GetPedProfileLowGain(i), 1.0, 1.0);
471     ((TProfile2D*)fPedestalHighGainRatio[i])->Divide(GetPedProfileHighGain(i), fReference->GetPedProfileHighGain(i), 1.0, 1.0);
472     ((TProfile2D*)fPeakMinusPedLowGainRatio[i])->Divide(GetPeakProfileLowGain(i), fReference->GetPeakProfileLowGain(i), 1.0, 1.0);
473     ((TProfile2D*)fPeakMinusPedHighGainRatio[i])->Divide(GetPeakProfileHighGain(i), fReference->GetPeakProfileHighGain(i), 1.0, 1.0);
474   
475     //For computing the difference, we cannot simply do TProfile2D->Add(), because that subtracts the sum of all entries,
476     //which means that the mean of the new profile will not be the difference of the means. So do it by hand:
477     for (int j = 0; j <= fColumns; j++) {
478       for (int k = 0; k <= fRows; k++) {
479         int bin = ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->GetBin(j+1, k+1);//Note that we assume here that all histos have the same structure...
480         double diff = fReference->GetPeakProfileHighGain(i)->GetBinContent(bin) - GetPeakProfileHighGain(i)->GetBinContent(bin);
481         ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->SetBinContent(j+1, k+1, diff);
482         ((TProfile2D*)fPeakMinusPedHighGainDiff[i])->SetBinEntries(bin, 1);
483
484         diff = fReference->GetPeakProfileLowGain(i)->GetBinContent(bin) - GetPeakProfileLowGain(i)->GetBinContent(bin);
485         ((TProfile2D*)fPeakMinusPedLowGainDiff[i])->SetBinContent(j+1, k+1, diff);
486         ((TProfile2D*)fPeakMinusPedLowGainDiff[i])->SetBinEntries(bin, 1);
487     
488         diff = fReference->GetPedProfileHighGain(i)->GetBinContent(bin) - GetPedProfileHighGain(i)->GetBinContent(bin);
489         ((TProfile2D*)fPedestalHighGainDiff[i])->SetBinContent(j+1, k+1, diff);
490         ((TProfile2D*)fPedestalHighGainDiff[i])->SetBinEntries(bin, 1);
491
492         diff = fReference->GetPedProfileLowGain(i)->GetBinContent(bin) - GetPedProfileLowGain(i)->GetBinContent(bin);
493         ((TProfile2D*)fPedestalLowGainDiff[i])->SetBinContent(j+1, k+1, diff);
494         ((TProfile2D*)fPedestalLowGainDiff[i])->SetBinEntries(bin, 1);
495        
496       } // rows
497     } // columns
498     
499   } // modules
500  
501 }
502
503 //_____________________________________________________________________
504 void AliCaloCalibPedestal::ComputeDeadTowers(int threshold, const char * deadMapFile)
505 {
506   //Computes the number of dead towers etc etc into memory, after this you can call the GetDead... -functions
507   int countTot = 0;
508   int countNew = 0;
509   int countRes = 0;
510   ofstream * fout = 0;
511   ofstream * diff = 0;
512   char name[512];//Quite a long temp buffer, just in case the filename includes a path
513   
514   if (deadMapFile) {
515     snprintf(name, 512, "%s.txt", deadMapFile);
516     fout = new ofstream(name);
517     snprintf(name, 512, "%sdiff.txt", deadMapFile);
518     diff = new ofstream(name);
519     if (!fout->is_open()) {
520       delete fout;
521       fout = 0;//Set the pointer to empty if the file was not opened
522     }
523     if (!diff->is_open()) {
524       delete diff;
525       fout = 0;//Set the pointer to empty if the file was not opened
526     }
527   }
528  
529   for (int i = 0; i < fModules; i++) {
530     if (GetPeakProfileHighGain(i)->GetEntries() > 0) { //don't care about empty histos
531       for (int j = 1; j <= fColumns; j++) {
532         for (int k = 1; k <= fRows; k++) {
533
534           if (GetPeakProfileHighGain(i)->GetBinContent(j, k) < threshold) {//It's dead
535             countTot++;//One more dead total
536             if (fout) {
537               (*fout) << i << " " 
538                       << (fRows - k) << " " 
539                       << (j-1) << " " 
540                       << "1" << " " 
541                       << "0" << endl;//Write the status to the deadmap file, if the file is open.
542             }
543             
544             if (fReference && fReference->GetPeakProfileHighGain(i)->GetBinContent(j, k) >= threshold) {
545               ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kRecentlyDeceased); 
546               countNew++;//This tower wasn't dead before!
547               if (diff) {
548                 ( *diff) << i << " " 
549                          << (fRows - k) << " " 
550                          << (j - 1) << " " 
551                          << "1" << " " 
552                          << "0" << endl;//Write the status to the deadmap difference file, if the file is open.
553               }
554             } 
555             else {
556               ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kDead);//This has been dead before. Nothing new         
557             }
558           } 
559           else { //It's ALIVE!!
560             //Don't bother with writing the live ones.
561             //if (fout)
562             //  (*fout) << i << " " 
563             //     << (fRows - k) << " " 
564             //     << (j - 1) << " " 
565             //     << "1" << " " 
566             //     << "1" << endl;//Write the status to the deadmap file, if the file is open.
567             if (fReference && fReference->GetPeakProfileHighGain(i)->GetBinContent(j, k) < threshold) {
568               ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kResurrected);
569               countRes++; //This tower was dead before => it's a miracle! :P
570               if (diff) {
571                 (*diff) << i << " " 
572                         << (fRows - k) << " " 
573                         << (j - 1) << " " 
574                         << "1" << " " 
575                         << "1" << endl;//Write the status to the deadmap difference file, if the file is open.
576               }
577             } 
578             else {
579               ((TH2D*)fDeadMap[i])->SetBinContent(j, k, kAlive);
580             }
581           }
582             
583         }//end for k/rows
584       }//end for j/columns
585     }//end if GetEntries >= 0
586   
587   }//end for modules
588  
589  if (fout) {
590    fout->close();
591    delete fout;
592  }
593  
594  fDeadTowers = countTot;
595  fNewDeadTowers = countNew;
596  fResurrectedTowers = countRes;
597 }
598