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