]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hshuttle.C
Adding comments (Laurent)
[u/mrichter/AliRoot.git] / HMPID / Hshuttle.C
1 void Hshuttle(Int_t runTime=1500)
2 {// this macro is to simulate the functionality of SHUTTLE.
3   gSystem->Load("libTestShuttle.so");
4   AliTestShuttle::SetMainCDB(TString("local://$HOME/CDB"));
5   
6   TMap *pDcsMap = new TMap;       pDcsMap->SetOwner(1);          //DCS archive map
7   
8   AliTestShuttle* pShuttle = new AliTestShuttle(0,0,1000000);   
9   SimPed();   for(Int_t ldc=1;ldc<=4;ldc++) pShuttle->AddInputFile(AliTestShuttle::kDAQ,"HMP","pedestals",Form("LDC%i",ldc),Form("HmpidPeds%i.tgz",ldc));
10   SimMap(pDcsMap,runTime); pShuttle->SetDCSInput(pDcsMap);                                    //DCS map
11   
12   AliPreprocessor* pp = new AliHMPIDPreprocessor(pShuttle); pShuttle->Process();  delete pp;  //here goes preprocessor 
13
14   DrawInput(pDcsMap); DrawOutput();
15 }//Hshuttle()
16 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
17 void SimPed()
18 {
19   ofstream out;
20   for(Int_t ddl=0;ddl<=13;ddl++){
21     out.open(Form("HmpidPedDdl%02i.txt",ddl));
22     out << 3 <<endl;
23     for(Int_t row=1;row<=24;row++)
24       for(Int_t dil=1;dil<=10;dil++)
25         for(Int_t adr=0;adr<=47;adr++){
26           Float_t mean  = 150+200*gRandom->Rndm();
27           Float_t sigma = ddl+gRandom->Gaus();
28           Int_t inhard=((Int_t(mean))<<9)+Int_t(mean+3*sigma);
29           out << Form("%2i %2i %2i %5.2f %5.2f %x\n",row,dil,adr,mean,sigma,inhard);
30         }
31
32     Printf("file ped %02i created",ddl);
33     out.close();
34   }
35   gSystem->Exec("tar czf HmpidPeds1.tgz HmpidPedDdl00.txt HmpidPedDdl01.txt HmpidPedDdl02.txt HmpidPedDdl03.txt");
36   gSystem->Exec("tar czf HmpidPeds2.tgz HmpidPedDdl04.txt HmpidPedDdl05.txt HmpidPedDdl06.txt HmpidPedDdl07.txt");
37   gSystem->Exec("tar czf HmpidPeds3.tgz HmpidPedDdl08.txt HmpidPedDdl09.txt HmpidPedDdl10.txt HmpidPedDdl11.txt");
38   gSystem->Exec("tar czf HmpidPeds4.tgz HmpidPedDdl12.txt HmpidPedDdl13.txt");
39   gSystem->Exec("rm -rf ped*.txt");
40 }
41 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
42 void SimMap(TMap *pDcsMap,Int_t runTime=1500)
43 {
44   Int_t stepTime=100; //time interval between mesuraments
45   Int_t startTime=0;
46   
47   
48   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
49     TObjArray *pP=new TObjArray;  pP->SetOwner(1);
50     TObjArray *pHV=new TObjArray; pHV->SetOwner(1); 
51     for(Int_t time=0;time<runTime;time+=stepTime)  pP->Add(new AliDCSValue((Float_t)1005.0 ,time));   //sample CH4 pressure [mBar]
52                                                    pHV->Add(new AliDCSValue((Float_t)2010.0,time));   //sample chamber HV [V]
53     pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_MP%i/HMP_MP%i_GAS/HMP_MP%i_GAS_PMWC.actual.value"           ,iCh,iCh,iCh)),pP); 
54     pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_MP%i/HMP_MP%i_PW/HMP_MP%i_SEC0/HMP_MP%i_SEC0_HV.actual.vMon",iCh,iCh,iCh)),pHV); 
55         
56     for(Int_t iRad=0;iRad<3;iRad++){//radiators loop
57       TObjArray *pT1=new TObjArray; pT1->SetOwner(1); 
58       TObjArray *pT2=new TObjArray; pT2->SetOwner(1); 
59       for (Int_t time=0;time<runTime;time+=stepTime)  pT1->Add(new AliDCSValue(13,time));  //sample inlet temperature    Nmean=1.292 @ 13 degrees
60       for (Int_t time=0;time<runTime;time+=stepTime)  pT2->Add(new AliDCSValue(14,time));  //sample outlet temperature
61       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iIn_Temp",iCh,iCh,iRad)) ,pT1); 
62       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iOut_Temp",iCh,iCh,iRad)),pT2);
63     }//radiators loop    
64   }//chambers loop
65 }//SimMap()
66 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
67 void DrawInput(TMap *pDcsMap)
68 {
69   TCanvas *c=new TCanvas("cc","Input data",600,600);    c->Divide(3,3);
70   
71   AliDCSValue *pVal; Int_t cnt;
72   
73   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
74     if(iCh==6) c->cd(1);  if(iCh==5) c->cd(2);                          
75     if(iCh==4) c->cd(4);  if(iCh==3) c->cd(5);  if(iCh==2) c->cd(6);
76                           if(iCh==1) c->cd(8);  if(iCh==0) c->cd(9); 
77                           
78     TObjArray *pHV=(TObjArray*)pDcsMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_PW/HMP_MP%i_SEC0/HMP_MP%i_SEC0_HV.actual.vMon",iCh,iCh,iCh,iCh)); //HV
79     TObjArray *pP =(TObjArray*)pDcsMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_GAS/HMP_MP%i_GAS_PMWC.actual.value",iCh,iCh,iCh)); //P
80     TGraph *pGrHV=new TGraph; pGrHV->SetMarkerStyle(5); TIter nextHV(pHV);
81     TGraph *pGrP =new TGraph; pGrP ->SetMarkerStyle(5); TIter nextP (pP );
82     
83     for(Int_t iRad=0;iRad<3;iRad++){
84       TObjArray *pT1=(TObjArray*)pDcsMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iIn_Temp",iCh,iCh,iRad));  TIter nextT1(pT1);
85       TObjArray *pT2=(TObjArray*)pDcsMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iOut_Temp",iCh,iCh,iRad)); TIter nextT2(pT2);
86       TGraph *pGrT1=new TGraph; pGrT1->SetMarkerStyle(5); 
87       TGraph *pGrT2=new TGraph; pGrT2->SetMarkerStyle(5); 
88       cnt=0; while((pVal=(AliDCSValue*)nextT1())) pGrT1->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());
89       cnt=0; while((pVal=(AliDCSValue*)nextT2())) pGrT2->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());
90       pGrT1->Draw("AP");  pGrT2->Draw("same");
91     }//radiators loop
92   }//chambers loop  
93 }//DrawInput()
94 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
95 void DrawOutput()
96 {
97   AliCDBManager::Instance()->SetDefaultStorage("local://$HOME/CDB"); AliCDBManager::Instance()->SetRun(0);
98   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre");
99   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean");
100   AliCDBEntry *pDaqSigEnt=AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");
101   
102   if(!pQthreEnt || !pNmeanEnt || !pDaqSigEnt) return;
103   
104   TObjArray *pNmean =(TObjArray*)pNmeanEnt ->GetObject(); 
105   TObjArray *pQthre =(TObjArray*)pQthreEnt ->GetObject(); 
106   TObjArray *pDaqSig=(TObjArray*)pDaqSigEnt->GetObject();
107    
108   TF1 *pRad0,*pRad1,*pRad2;  
109   TCanvas *c2=new TCanvas("c2","Output"); c2->Divide(3,3);
110   
111   TH1F *pSig[7];
112   
113   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
114     TMatrix *pM=(TMatrix*)pDaqSig->At(iCh);
115     
116     pSig[iCh]=new TH1F(Form("sig%i",iCh),"Sigma;[QDC]",100,-5,20); //pSig[iCh]->SetLineColor(iCh+kRed);
117     for(Int_t padx=0;padx<160;padx++) for(Int_t pady=0;pady<144;pady++) pSig[iCh]->Fill((*pM)(padx,pady));
118     
119     c2->cd(7);    if(iCh==0) pSig[iCh]->Draw(); else pSig[iCh]->Draw("same");
120     
121     if(iCh==6) c2->cd(1);  if(iCh==5) c2->cd(2);                          
122     if(iCh==4) c2->cd(4);  if(iCh==3) c2->cd(5);  if(iCh==2) c2->cd(6);
123                            if(iCh==1) c2->cd(8);  if(iCh==0) c2->cd(9); 
124                           
125     TF1 *pRad0=(TF1*)pNmean->At(iCh*3+0); pRad0->Draw();  pRad0->GetXaxis()->SetTimeDisplay(kTRUE); pRad0->GetYaxis()->SetRangeUser(1.28,1.3);
126     TF1 *pRad1=(TF1*)pNmean->At(iCh*3+1); pRad1->Draw("same");
127     TF1 *pRad2=(TF1*)pNmean->At(iCh*3+2); pRad2->Draw("same");
128   }//chambers loop    
129 }//DrawOutput()
130 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++