]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackerHV.cxx
fix finding of pad neighbours; remove methods to write them in OCDB
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackerHV.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
16 // $Id$
17
18 #include "AliMUONTrackerHV.h"
19
20 #include <algorithm>
21 #include <map>
22 #include <set>
23
24 #include "AliCDBManager.h"
25 #include "AliCDBEntry.h"
26 #include "AliDCSValue.h"
27 #include "AliGRPObject.h"
28 #include "AliMpArrayI.h"
29 #include "AliMpConstants.h"
30 #include "AliMpDCSNamer.h"
31 #include "AliMpDEStore.h"
32 #include "AliMpDetElement.h"
33 #include "AliMUON2DMap.h"
34 #include "AliMUONCalibParamND.h"
35 #include "AliMUONCalibrationData.h"
36 #include "AliMUONCDB.h"
37 #include "AliMUONPainterDataRegistry.h"
38 #include "AliMUONTrackerData.h"
39 #include "AliMUONTrackerDataWrapper.h"
40 #include "AliLog.h"
41
42 #include "TCanvas.h"
43 #include "TGraph.h"
44 #include "TH2.h"
45 #include "TLine.h"
46 #include "TMap.h"
47 #include "TMultiGraph.h"
48 #include "TObjArray.h"
49 #include "TObjString.h"
50 #include "TStyle.h"
51 #include "Riostream.h"
52
53 //
54 // Class to inspect the MUON TRACKER HV values
55 //
56 // With this class you can :
57 //
58 // a) get a list of trips (method ReportTrips)
59 // b) print the values for some (or all) HV channels (method Print)
60 // c) plot the values for some (or all) HV channels (method Plot)
61 // d) get a list of HV channels that are "OFF" (methods Scan and HVoff)
62 //
63 // Note that in this class, all the output (either text or canvas) or the
64 // channel *names* used are the same as in the DCS UI at Pt2
65 // Specifically the chamber ids start at 1, the slat numbers at 1 and
66 // the quad and sect number at 1 also. And not at zero like for the
67 // DCS *aliases*. On the contraty, the internal map, coming from the OCDB,
68 // only contains aliases, not names. Confusing ? It is.
69 //
70
71 ClassImp(AliMUONTrackerHV)
72
73 //______________________________________________________________________________
74 AliMUONTrackerHV::AliMUONTrackerHV(const char* runlist, const char* ocdbPath)
75 : TObject(), fRunList(), fOCDBPath(ocdbPath), fDCSNamer(0x0)
76 {
77   // ctor from a runlist (txt file)
78   SetRunList(runlist);
79 }
80
81 //______________________________________________________________________________
82 AliMUONTrackerHV::AliMUONTrackerHV(Int_t runNumber, const char* ocdbPath)
83 : TObject(), fRunList(), fOCDBPath(ocdbPath), fDCSNamer(0x0)
84 {
85   // ctor for a single run
86   SetRunList(runNumber);
87 }
88
89 //______________________________________________________________________________
90 AliMUONTrackerHV::~AliMUONTrackerHV()
91 {
92   // dtor
93   delete fDCSNamer;
94 }
95
96 //____________________________________________________________________________
97 TMultiGraph* AliMUONTrackerHV::CombineMulti(TObjArray& graphs)
98 {
99   // combine multigraphs
100   
101   TMultiGraph* rv = new TMultiGraph;
102   
103   TIter next(&graphs);
104   TMultiGraph* mg;
105   TMultiGraph* ref = static_cast<TMultiGraph*>(next());
106
107   Int_t dref = ref->GetListOfGraphs()->GetEntries();
108
109   while ( ( mg = static_cast<TMultiGraph*>(next())) )
110   {
111     TList* list = mg->GetListOfGraphs();
112     Int_t d1 = list->GetEntries();
113     
114     if (  d1 != dref )
115     {
116       AliError(Form("%d vs %d",d1,dref));
117       return 0x0;
118     }
119   }
120   
121   for ( Int_t i = 0; i < dref; ++i )
122   {    
123     TObjArray graph;
124     next.Reset();
125     while ( ( mg = static_cast<TMultiGraph*>(next())) )
126     {
127       graph.Add(mg->GetListOfGraphs()->At(i));
128       TGraph* g = Combine(graph);
129       rv->Add(g);
130     }
131   }
132   return rv;
133 }
134
135 //____________________________________________________________________________
136 TGraph* AliMUONTrackerHV::Combine(TObjArray& graphs)
137 {
138   // make one graph out of several
139   // x axis is supposed to be time and will end up ordered in the
140   // returned graph
141   
142   std::map<int, std::vector<double> > values;
143   std::map<int, std::vector<double> >::const_iterator it;
144   
145   TIter next(&graphs);
146   TGraph* g;
147   
148   while ( ( g = static_cast<TGraph*>(next())) )
149   {
150     for ( Int_t i = 0; i < g->GetN(); ++i )
151     {
152       std::vector<double> pair;
153       
154       pair.push_back(g->GetX()[i]);
155       pair.push_back(g->GetY()[i]);
156       
157       values.insert( std::make_pair(g->GetX()[i],pair));
158     }
159   }
160   
161   TGraph* rv(0x0);
162   
163   if ( values.size() )
164   {
165     std::vector<double> vx;
166     std::vector<double> vy;
167     
168     for ( it = values.begin(); it != values.end(); ++it )
169     {
170       const std::vector<double>& q = it->second;
171       
172       vx.push_back(q[0]);
173       vy.push_back(q[1]);
174     }
175     
176     rv = new TGraph(values.size(),&vx[0],&vy[0]);
177     rv->GetXaxis()->SetNoExponent();
178     
179     g = static_cast<TGraph*>(graphs.At(0));
180     
181     rv->SetName(g->GetName());
182     rv->SetTitle(g->GetTitle());
183   }
184   
185   return rv;
186 }
187
188 //______________________________________________________________________________
189 void AliMUONTrackerHV::ReadIntegers(const char* filename, std::vector<int>& integers)
190 {
191   /// Read integers from filename, where integers are either
192   /// separated by "," or by return carriage
193   std::ifstream in(gSystem->ExpandPathName(filename));
194   int i;
195   
196   std::set<int> runset;
197   
198   char line[10000];
199   
200   in.getline(line,10000,'\n');
201   
202   TString sline(line);
203   
204   if (sline.Contains(","))
205   {
206     TObjArray* a = sline.Tokenize(",");
207     TIter next(a);
208     TObjString* s;
209     while ( ( s = static_cast<TObjString*>(next()) ) )
210     {
211       runset.insert(s->String().Atoi());
212     }
213     delete a;
214   }
215   else
216   {
217     runset.insert(sline.Atoi());
218     
219     while ( in >> i )
220     {
221       runset.insert(i);
222     }
223   }
224   
225   for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
226   {
227     integers.push_back((*it));
228   }
229   
230   std::sort(integers.begin(),integers.end());
231 }
232
233 //______________________________________________________________________________
234 AliMpDCSNamer*
235 AliMUONTrackerHV::DCSNamer() const
236 {
237   // return the dcs namer
238   if (!fDCSNamer)
239   {
240     if (!AliMpDEStore::Instance(false))
241     {
242       AliMUONCDB::LoadMapping();
243     }
244     fDCSNamer = new AliMpDCSNamer("TRACKER");
245   }
246   return fDCSNamer;
247 }
248
249 //______________________________________________________________________________
250 void AliMUONTrackerHV::SetRunList(Int_t runNumber)
251 {
252   // Make the runlist be a single run
253   fRunList.clear();
254   fRunList.push_back(runNumber);
255 }
256
257 //______________________________________________________________________________
258 void
259 AliMUONTrackerHV::SetRunList(const char* runlist)
260 {
261   // Read the runlist from an ASCII file or a comma separated list
262   // or a space separated list
263   
264   fRunList.clear();
265   
266   if ( TString(runlist).Contains(",") || TString(runlist).Contains(" ") )
267   {
268     TObjArray* runs = 0x0;
269     if ( TString(runlist).Contains(",") )
270     {
271       runs = TString(runlist).Tokenize(",");
272     }
273     else
274     {
275       runs = TString(runlist).Tokenize(" ");
276     }
277     TIter next(runs);
278     TObjString* s;
279     std::set<int> runset;
280     
281     while ( ( s = static_cast<TObjString*>(next()) ) )
282     {
283       runset.insert(s->String().Atoi());
284     }
285     
286     for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
287     {
288       fRunList.push_back((*it));
289     }
290     
291     std::sort(fRunList.begin(),fRunList.end());
292     
293     delete runs;
294   }
295   else
296   {
297     ReadIntegers(runlist,fRunList);
298   }
299 }
300
301
302 //______________________________________________________________________________
303 TGraph*
304 AliMUONTrackerHV::GraphValues(TMap* m, const char* dcsname)
305 {
306   // make a graph of HV channels' voltage values for a given dcs name (name, not
307   // alias)
308   
309   if ( TString(dcsname).Contains("sw") )
310   {
311     // do not graph switches
312     return 0x0;
313   }
314
315   
316   AliInfo(dcsname);
317   
318   TPair* p = static_cast<TPair*>(m->FindObject(DCSNamer()->DCSAliasFromName(dcsname).Data()));
319   
320   if (!p) return 0x0;
321   
322   TObjArray* a = static_cast<TObjArray*>(p->Value());
323   TIter n2(a);
324   AliDCSValue* val;
325   Int_t i(0);
326
327   TGraph* g = new TGraph(a->GetEntries());
328   while ( ( val = static_cast<AliDCSValue*>(n2()) ) )
329   {
330     g->SetPoint(i,val->GetTimeStamp(),val->GetFloat());
331     ++i;
332   }
333   g->SetName(dcsname);
334   return g;
335 }
336
337 //______________________________________________________________________________
338 void
339 AliMUONTrackerHV::Scan(Int_t verbose)
340 {
341   /// Retrieve HV values from OCDB for a given run list, and check whether
342   /// we have some issues with them...
343   /// If you pipe the results of this into a text file, you can then
344   /// feed it to the HVoff method for further investigations.
345   ///
346   
347   if ( fRunList.empty() )
348   {
349     AliError("No runs to process...");
350     return;    
351   }
352     
353   AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
354   
355   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
356   {
357     AliMUONCDB::CheckHV(fRunList[i],verbose);
358   }
359 }
360
361 //______________________________________________________________________________
362 void AliMUONTrackerHV::HVoff(const char* logfile, const char* outputBaseName)
363 {
364   /// Check the number of HV which have problem
365   /// the input is the output of e.g.
366   /// .L MUONTrackerHV.C+
367   /// ScanHV("lhc11de.list");>  lhc11de.log
368   ///
369   
370   gStyle->SetOptStat(0);
371   
372   char line[1024];
373   
374   std::ifstream in(logfile);
375   int run(-1),a,b,c,d,e,f,g,h,z,other;
376   std::map<int,std::string> results;
377   
378   std::string message;
379   const char* testProblem = "I-AliMUONCDB::CheckHV::CheckHV:      Problem at ";
380   
381   while ( in.getline(line,1023,'\n') )
382   {
383     TString sline(line);
384     if (sline.Contains("SUMMARY"))
385     {
386       AliInfo(line);
387       int r;
388       sscanf(line,"I-AliMUONCDB::CheckHV::CheckHV: RUN %09d HVchannel SUMMARY : # of cases A(%3d) B(%3d) C(%3d) D(%3d) E(%3d) F(%3d) G(%3d) H(%3d) Z(%3d) OTHER(%3d)",
389              &r,&a,&b,&c,&d,&e,&f,&g,&h,&z,&other);
390       if ( r != run )
391       {
392         if ( run == -1 )
393         {
394           AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
395           AliCDBManager::Instance()->SetRun(r);
396           AliMUONCDB::LoadMapping();
397         }
398         
399         if ( run > 0 )
400         {
401           results.insert(std::make_pair(run,message));
402           
403         }
404         message = "";
405         run = r;
406       }          
407     }
408     else if ( sline.Contains(testProblem) )
409     {
410       message += "|";
411       message += sline(strlen(testProblem),sline.Length()-1).Data();
412     }
413   }
414   
415   results.insert(std::make_pair(run,message));
416   
417   TH2* hvoff = new TH2I(outputBaseName,outputBaseName,1,0,1,1,0,1);
418   
419   std::map<int,std::string>::const_iterator it;
420   
421   for ( it = results.begin(); it != results.end(); ++it )
422   {
423     AliInfo(Form("%d -> %s",it->first,it->second.c_str()));
424     TObjArray* split = TString(it->second.c_str()).Tokenize("|");
425     TIter next(split);
426     TObjString* str;
427     while ( ( str = static_cast<TObjString*>(next()) ) )
428     {
429       TString s(str->String());
430       TObjArray* parts = s.Tokenize(":");
431       TString alias = (static_cast<TObjString*>(parts->At(0)))->String();
432       TString channel = DCSNamer()->DCSNameFromAlias(alias.Data());
433       channel += Form("(%4d)",DCSNamer()->DetElemIdFromDCSAlias(alias.Data()));
434       channel.ReplaceAll(".actual.vMon","");
435       hvoff->Fill(Form("%6d",it->first),channel.Data(),1.0);
436       delete parts;
437     }
438     delete split;
439   }
440   
441   hvoff->LabelsDeflate("x");
442   hvoff->LabelsDeflate("y");
443   hvoff->LabelsOption("x","<");
444   hvoff->LabelsOption("y","<");
445   
446   TCanvas* c1 = new TCanvas;
447   c1->SetLeftMargin(0.35);
448   hvoff->Draw("text");
449   c1->Print(Form("%s.pdf",outputBaseName));
450   TCanvas* c2 = new TCanvas;
451   TH1* hx = hvoff->ProjectionX("hvoffperrun");
452   hx->Draw();
453   c2->Print(Form("%s-perrun.pdf",outputBaseName));
454   TCanvas* c3 = new TCanvas;
455   c3->SetBottomMargin(0.55);
456   TH1* perchannel = hvoff->ProjectionY("hvoffperchannel");
457   perchannel->GetXaxis()->SetBit(TAxis::kLabelsVert);
458   perchannel->GetXaxis()->LabelsOption(">");
459   perchannel->Draw("texthist");
460   c3->Print(Form("%s-perchannel.pdf",outputBaseName));
461 }
462
463 //______________________________________________________________________________
464 void AliMUONTrackerHV::TimeAxis(TMultiGraph* g)
465 {
466   g->GetXaxis()->SetTimeDisplay(1);
467 //  g->GetXaxis()->SetTimeFormat("%d/%m %H:%M%F2010-12-31 24:00:00");
468   g->GetXaxis()->SetTimeFormat("%d/%m %H:%M");
469   g->GetXaxis()->SetTimeOffset(0,"gmt");
470   g->GetXaxis()->SetNdivisions(505);
471 }
472
473 //______________________________________________________________________________
474 TMultiGraph*
475 AliMUONTrackerHV::GraphHV(TMap* m, const char* dcsname)
476 {
477   // Make a graph of the values matching dcsname
478   TIter next(m);
479   TObjString* s;
480   
481   TMultiGraph* mg = new TMultiGraph;
482
483   while ( ( s = static_cast<TObjString*>(next()) ) )
484   {
485     TString name(DCSNamer()->DCSNameFromAlias(s->String()));
486     
487     if ( dcsname && !name.Contains(dcsname)) continue;
488     
489     TGraph* g = GraphValues(m,name);
490     
491     if ( g ) 
492     {
493       g->SetMarkerSize(1.5);
494       g->SetMarkerStyle(2);
495       g->SetLineStyle(2);
496       mg->Add(g,"lp");
497       g->SetTitle(name.Data());
498     }
499   }  
500
501   return mg;
502 }
503
504 //______________________________________________________________________________
505 void
506 AliMUONTrackerHV::Print(Option_t* dcsname) const
507 {
508   /// Print HV values for a given dcs name (or all if dcsname=0)
509   
510   AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
511   TList messages;
512   messages.SetOwner(kTRUE);
513   
514   for ( std::vector<int>::size_type iRun = 0; iRun < fRunList.size(); ++iRun )
515   {
516     Int_t runNumber = fRunList[iRun];
517     
518     AliInfo("---------------------");
519     AliInfo(Form("RUN %09d",runNumber));
520     
521     messages.Delete();
522     
523     AliCDBManager::Instance()->SetRun(runNumber);
524     
525     Bool_t patchValues(kFALSE);
526     Bool_t dryRun(kTRUE);
527     
528     TMap* m = AliMUONCalibrationData::CreateHV(runNumber,0x0,patchValues,&messages,dryRun);
529     
530     TIter next(m);
531     TObjString* s;
532     
533     while ( ( s = static_cast<TObjString*>(next()) ) )
534     {
535       TString name(DCSNamer()->DCSNameFromAlias(s->String()));
536       
537       if ( dcsname && !name.Contains(dcsname)) continue;
538
539       TPair* p = static_cast<TPair*>(m->FindObject(s->String()));
540       
541       if (!p) continue;
542       
543       TObjArray* a = static_cast<TObjArray*>(p->Value());
544       TIter n2(a);
545       AliDCSValue* val;
546       Int_t i(0);
547       
548       while ( ( val = static_cast<AliDCSValue*>(n2()) ) )
549       {
550         std::cout << Form("i=%5d ",i) << std::endl;
551         val->Print("");
552         ++i;
553       }
554     }
555   }
556 }
557
558 //______________________________________________________________________________
559 void
560 AliMUONTrackerHV::Plot(const char* dcsname, Bool_t withPatch, Bool_t plotIntermediate)
561 {
562   /// Show HV values for a given dcs name (or all if dcsname=0)
563   /// Each canvas for each run will go to a separate PDF file
564   
565   AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
566   TList messages;
567   messages.SetOwner(kTRUE);
568   TObjArray graphs;
569   
570   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
571   {
572     Int_t runNumber = fRunList[i];
573   
574     messages.Delete();
575     
576     AliCDBManager::Instance()->SetRun(runNumber);
577     
578     TMap* m = AliMUONCalibrationData::CreateHV(runNumber,0x0,withPatch,&messages,kTRUE);
579     
580     TMultiGraph* mg = GraphHV(m,dcsname);
581     
582     if ( !mg ) continue;
583     
584     graphs.Add(mg);
585     
586     TString cname(Form("MCH_HV_RUN%09d",runNumber));
587     
588     if ( strlen(dcsname) > 0 )
589     {
590       TString s(dcsname);
591       s.ReplaceAll("/","_");
592       cname += Form("_dcsname_%s",s.Data());
593     }
594     
595     AliCDBEntry* e = AliCDBManager::Instance()->Get("GRP/GRP/Data",runNumber);
596     
597     TLine* startRunLine(0);
598     TLine* endRunLine(0);
599     time_t start(0);
600     time_t end(0);
601     
602     if ( e )
603     {
604       AliGRPObject* grp = static_cast<AliGRPObject*>(e->GetObject());
605       if (grp)
606       {
607         start = grp->GetTimeStart();
608         end = grp->GetTimeEnd();
609       }
610     }
611     
612     if ( end )
613     {
614       TGraph* g = new TGraph(1);
615       g->SetPoint(0,end,0);
616       mg->Add(g,"");
617     }
618     
619     if ( plotIntermediate )
620     {
621       TCanvas* c = new TCanvas(cname.Data(),cname.Data());
622     
623       c->Draw();
624     
625       mg->SetTitle(cname.Data());
626     
627       mg->Draw("AL");
628     
629       TimeAxis(mg);
630     
631       if ( start )
632       {
633         startRunLine = new TLine(start,mg->GetYaxis()->GetXmin(),start,mg->GetYaxis()->GetXmax());
634         startRunLine->SetLineColor(2);
635         startRunLine->SetLineWidth(4);
636       }
637       if  ( end )
638       {
639         endRunLine = new TLine(end,mg->GetYaxis()->GetXmin(),end,mg->GetYaxis()->GetXmax());
640         endRunLine->SetLineColor(2);
641         endRunLine->SetLineWidth(4);
642       }
643     
644       if ( startRunLine ) startRunLine->Draw();
645       if ( endRunLine ) endRunLine->Draw();
646     
647       c->SaveAs(Form("%s.pdf",cname.Data()));
648     }
649   }
650   
651   new TCanvas;
652   
653   TMultiGraph* g = CombineMulti(graphs);
654
655   TIter next(g->GetListOfGraphs());
656   TGraph* gi;
657   
658   while ( ( gi = static_cast<TGraph*>(next())))
659   {
660     gi->SetMarkerStyle(kPlus);
661   }
662   g->Draw("alp");
663   TimeAxis(g);
664 }
665
666 //______________________________________________________________________________
667 void
668 AliMUONTrackerHV::ReportTrips(Bool_t includeLowOnes)
669 {
670   /// Report trips
671   /// if includeLowOnes is kTRUE we'll report also the trips which starts from non-operational voltage values
672   
673   AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
674   
675   TList messages;
676   messages.SetOwner(kTRUE);
677   TObjString* msg(0);
678
679   std::map<std::string,int> channels;
680
681   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
682   {
683     Int_t runNumber = fRunList[i];
684     
685     AliInfo("---------------------");
686     
687     Int_t ntrips(0);
688     
689     messages.Delete();
690     
691     AliCDBManager::Instance()->SetRun(runNumber);
692     
693     AliMUONCalibrationData::CreateHV(runNumber,0x0,kTRUE,&messages,kTRUE);
694     
695     if (!AliMpDEStore::Instance(false))
696     {
697       AliMUONCDB::LoadMapping();
698     }
699     
700     TIter next(&messages);
701
702     while ( ( msg = static_cast<TObjString*>(next())) )
703     {
704       if ( msg->String().Contains("TRIP") && ( includeLowOnes || !msg->String().Contains("LOWTRIP") ) )
705       {
706         ++ntrips;
707       }
708     }
709
710     AliInfo(Form("RUN %09d - %d trip%c",runNumber,ntrips,(ntrips>1 ? 's':' ')));
711     
712     next.Reset();
713     std::map<int,std::string> report;
714     
715     while ( ( msg = static_cast<TObjString*>(next())) )
716     {
717       if ( msg->String().Contains("TRIP") )
718       {
719         TObjArray* parts = msg->String().Tokenize(" ");
720         TString channelName(static_cast<TObjString*>(parts->At(0))->String());
721         
722         for ( Int_t ip = 0; ip <= parts->GetLast(); ++ip)
723         {
724           TString p(static_cast<TObjString*>(parts->At(ip))->String());
725           
726           if ( p.Contains("TRIP") )
727           {
728             if ( includeLowOnes || !p.Contains("LOWTRIP") )
729             {
730               TString ts(static_cast<TObjString*>(parts->At(ip+2))->String());
731           
732               ip += 3;
733           
734               Int_t index = ts.Index("TS:");
735           
736               UInt_t timeStamp = TString(ts(index+strlen("TS:"),ts.Length()-index)).Atoi();
737           
738               TString tmp(msg->String());
739               tmp.ReplaceAll(channelName.Data(),DCSNamer()->DCSNameFromAlias(channelName.Data()));
740               report[timeStamp] = tmp.Data();
741               channels[channelName.Data()]++;
742             }
743           }
744         }
745         delete parts;
746       }
747     }
748
749     for ( std::map<int,std::string>::const_iterator it = report.begin(); it != report.end(); ++it )
750     {
751       AliInfo(Form("%s %s",TTimeStamp(it->first).AsString("s"),it->second.c_str()));
752     }
753   }
754   
755   AliInfo("--------------------------------------------------------------------");
756   AliInfo("++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
757   
758   int totalTrips(0);
759   AliMUON2DMap tripMap(kTRUE);
760   Int_t nofChannels(AliMpConstants::ManuNofChannels());
761
762   for ( std::map<std::string,int>::const_iterator it = channels.begin(); it != channels.end(); ++it )
763   {
764     AliInfo(Form("%40s %3d",DCSNamer()->DCSNameFromAlias(it->first.c_str()).Data(),it->second));
765     totalTrips += it->second;
766     
767     Int_t detElemId = DCSNamer()->DetElemIdFromDCSAlias(it->first.c_str());
768     
769     AliMpDetElement* de = AliMpDEStore::Instance()->GetDetElement(detElemId);
770     
771     // build the list of manuIds for this channel
772     AliMpArrayI manuArray;
773     
774     manuArray.SetSize(300);
775     
776     Int_t index = DCSNamer()->DCSIndexFromDCSAlias(it->first.c_str());
777     Int_t firstIndex(index);
778     Int_t lastIndex(index);
779     
780     if ( index < 0 )
781     {
782       // it's a slat, must loop over PCBs
783       firstIndex = 0;
784       lastIndex = DCSNamer()->NumberOfPCBs(detElemId)-1;
785     }
786     
787     for ( int i = firstIndex; i <= lastIndex ; ++i )
788     {
789       const AliMpArrayI* ma = de->ManusForHV(i);
790       if (!ma)
791       {
792         AliError(Form("Could not get ma for de %d index %d",detElemId,i));
793         continue;
794       }
795       for ( int j = 0; j < ma->GetSize(); ++j )
796       {
797         manuArray.Add(ma->GetValue(j),kFALSE);
798       }
799     }
800     
801     for ( Int_t iManu = 0; iManu < manuArray.GetSize(); ++iManu )
802     {
803       Int_t manuId = manuArray.GetValue(iManu);
804       
805       AliMUONVCalibParam* tripRate = new AliMUONCalibParamND(1,nofChannels,detElemId,manuId,0);
806       
807       tripMap.Add(tripRate);
808       
809       for ( Int_t j = 0 ; j < nofChannels; ++j )
810       {
811         tripRate->SetValueAsDouble(j,0,it->second*1.0);
812       }
813     }
814   }
815
816   AliInfo("++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
817   AliInfo(Form("Total of %3d trips for %4ld runs",totalTrips,fRunList.size()));
818   
819   AliMUONTrackerData* data = new AliMUONTrackerData("tripcount","Number of trips",1);
820   data->Add(tripMap);
821   data->SetDimensionName(0,"ntrips");
822
823   AliMUONVTrackerDataMaker* dw = new AliMUONTrackerDataWrapper(data);
824   
825   AliMUONPainterDataRegistry::Instance()->Register(dw);
826
827 }
828