]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Comparing lists of noisy SSD chips (Panos)
authormasera <masera@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 21 Nov 2009 18:21:23 +0000 (18:21 +0000)
committermasera <masera@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 21 Nov 2009 18:21:23 +0000 (18:21 +0000)
ITS/readSSDCommonMode.C
ITS/readSSDOCDBEntry.C

index 3f5c6e75b5a637d2529e0e3519f3c87a875b8d60..d28a217c94c9aabd8a9fe0ed4d3ed5a4313e8665 100644 (file)
@@ -38,6 +38,9 @@ static const Int_t fgkNumberOfPSideStrips = 768; //number of P-side strips
 TList *initCM();
 void makeCM(const char* filename, Int_t nEvents, TList *list);
 void checkCM(const char* filename);
+void compareChipLists(TString inputFile1, 
+                     TString inputFile2, 
+                     TString outputTxt);
 
 //__________________________________________________________//
 void readSSDCommonMode(const char* filename = "raw.root",
@@ -434,3 +437,54 @@ void checkCM(const char* filename) {
 
   fInput->Close();
 }
+
+//__________________________________________________________//
+void compareChipLists(TString inputFile1, 
+                     TString inputFile2, 
+                     TString outputTxt){
+  //Compare two different lists of chips and put those present in both files
+  //in a new list
+  //open input lists
+  TFile *fInput1 = TFile::Open(inputFile1.Data());
+  TList *listOfSuspiciousChips1 = new TList();
+  listOfSuspiciousChips1=fInput1->GetListOfKeys();
+  TFile *fInput2 = TFile::Open(inputFile2.Data());
+  TList *listOfSuspiciousChips2 = new TList();
+  listOfSuspiciousChips2=fInput2->GetListOfKeys();
+       
+  Int_t Nentries1 = 0, Nentries2 = 0, k = 0;
+  Nentries1 = listOfSuspiciousChips1->GetEntries();
+  Nentries2 = listOfSuspiciousChips2->GetEntries();
+       
+  //create new list
+  TList *listOfRecurrentChips = new TList();
+  TString Name1;
+       
+  for(Int_t i=0; i<Nentries1; i++){
+    TH1F *h1 = dynamic_cast <TH1F*>(fInput1->Get(listOfSuspiciousChips1->At(i)->GetName()));
+    for(Int_t j=0; j<Nentries2; j++){
+      TH1F *h2 = dynamic_cast <TH1F*>(fInput2->Get(listOfSuspiciousChips2->At(j)->GetName()));
+      Name1=h1->GetName();
+      if (!Name1.CompareTo(h2->GetName())) {
+       cout << Name1.CompareTo(h2->GetName()) << 
+         "   " << h2->GetName() << endl;
+       k++;
+       listOfRecurrentChips->Add(h1);
+      }
+    }//second file loop
+  }//first file loop
+  
+  Printf("%i Faulty chips in the first file", Nentries1);
+  Printf("%i Faulty chips in the second file", Nentries2);
+  Printf("%i Recurrent Faulty chips", k+1);    
+  
+  TString outputFile = "SSD.RecurrentFaultyChips."; outputFile += outputTxt;
+  outputFile += ".root";
+  
+  TFile *fOutput = new TFile(outputFile.Data(),"recreate");
+  listOfRecurrentChips->Write();
+  fOutput->Close();
+  
+  fInput1->Close();
+  fInput2->Close();
+}
index 3f6db98c3172d25b72e1b8b599e73e29297e541e..6e65dd9696f7366f909eb8a9f842aecb09d1f5be 100644 (file)
@@ -466,6 +466,13 @@ void BadChannelMap(AliCDBManager * man) {
   cBadChannel->cd(4)->SetGridx(); cBadChannel->cd(4)->SetGridy();
   cBadChannel->cd(4); fHistNSideBadChannelMapLayer6->Draw("colz");
   cBadChannel->SaveAs("Run-BadChannels.gif");
+
+  TFile *fOutput = new TFile("badChannelsSSD.root","recreate");
+  fHistPSideBadChannelMapLayer5->Write();
+  fHistNSideBadChannelMapLayer5->Write();
+  fHistPSideBadChannelMapLayer6->Write();
+  fHistNSideBadChannelMapLayer6->Write();
+  fOutput->Close();
 }
 
 //_____________________________________________________________________//