]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hshuttle.C
Warning fixed
[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   AliTestShuttle::SetMainCDB(TString("local://$ALICE_ROOT/OCDB"));
7   
8   TMap *pDcsMap = new TMap;       pDcsMap->SetOwner(1);          //DCS archive map
9   
10   AliTestShuttle* pShuttle = new AliTestShuttle(0,0,1000000);   
11   //pShuttle->SetInputRunType("PHYSICS");
12   pShuttle->SetInputRunType("CALIBRATION");
13   //pShuttle->SetInputRunType("STANDALONE");
14   SimPed();   
15   for(Int_t ldc=51;ldc<=52;ldc++) 
16   {
17    if(ldc==51) {for(Int_t iddl=0;iddl<=7;iddl++)  pShuttle->AddInputFile(AliTestShuttle::kDAQ,"HMP",Form("HmpidPedDdl%02i.txt",iddl),Form("LDC%i",ldc),Form("HmpidPedDdl%02i.txt",iddl));}
18    if(ldc==52) {for(Int_t iddl=8;iddl<=13;iddl++) pShuttle->AddInputFile(AliTestShuttle::kDAQ,"HMP",Form("HmpidPedDdl%02i.txt",iddl),Form("LDC%i",ldc),Form("HmpidPedDdl%02i.txt",iddl));}
19   }
20   SimMap(pDcsMap,runTime); pShuttle->SetDCSInput(pDcsMap);                                    //DCS map
21   SimNoiseMap();           pShuttle->AddInputFile(AliTestShuttle::kDAQ,"HMP", "HmpPhysicsDaNoiseMap.root", "MON", "./HmpPhysicsDaNoiseMap.root");                                                                          //Simulate Noise Map
22   
23     
24   AliPreprocessor* pp = new AliHMPIDPreprocessor(pShuttle); pShuttle->Process();              //here goes preprocessor  
25 //  delete pp;   
26
27   DrawInput(pDcsMap); DrawOutput();
28   
29   gSystem->Exec("rm -rf HmpidPedDdl*.txt");                 //delete simulated files
30   gSystem->Exec("rm -rf ./HmpPhysicsDaNoiseMap.root");      //delete simulated file
31   
32 }//Hshuttle()
33 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 void SimNoiseMap()
35 {
36   //
37   // Simulate the noise map that will be extracted by the HMP Physics DA from PHYSICS and STANDALONE runs
38   //
39   TRandom rnd;
40   TH2F *hHmpNoiseMaps=new TH2F("hHmpNoiseMaps","Noise Maps Ch: 0-7;Ch 0-7: pad X;Ch0, Ch1, Ch2, Ch3, Ch4, Ch5, Ch6 pad Y ;Noise level (\%)",160,0,160,144*7,0,144*7); //In y we have all 7 chambers
41   Int_t events=2000;
42   for(Int_t nev=0;nev<events;nev++) 
43     for(Int_t npads=0;npads<350;npads++)  //average 5- pads times 7 chambers
44       hHmpNoiseMaps->Fill(rnd.Integer(160),rnd.Integer(7)*144+rnd.Integer(144));
45   
46   if(events!=0) hHmpNoiseMaps->Scale(1.0/events);
47   TFile *fout = new TFile("HmpPhysicsDaNoiseMap.root","recreate");
48   hHmpNoiseMaps->Write();
49   fout->Close();
50   Printf("==> HMP Noise Map simulation has finished!");
51 }
52 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
53 void SimPed()
54 {
55   Int_t iDDLmin=0,iDDLmax=13;
56   Int_t nSigmas = 3;                                            // value stored in the ddl files of pedestals
57   ofstream out;
58   for(Int_t ddl=iDDLmin;ddl<=iDDLmax;ddl++){
59     out.open(Form("HmpidPedDdl%02i.txt",ddl));
60     out << Form("%8s %2d\n","RunNumber",      999);             //read run number
61     out << Form("%8s %2d\n","LdcId" ,         999);             //read LDC Id
62     out << Form("%8s %2d\n","TimeStamp",      999);             //read time stamp
63     out << Form("%8s %2d\n","TotNumEvt",      999);             //read number of total events processed
64     out << Form("%8s %2d\n","TotDDLEvt",      999);             //read number of bad events for DDL # nDDL processed
65     out << Form("%8s %2d\n","NumBadEvt",      999);             //read number of bad events for DDL # nDDL processed
66     out << Form("%8s %2f\n","NBadE(%)",       999.9);           //read number of bad events (in %) for DDL # nDDL processed
67     out << Form("%8s %2.2d\n","SigCut",      nSigmas);          //# of sigma cuts
68     for(Int_t row=1;row<=24;row++)
69       for(Int_t dil=1;dil<=10;dil++)
70         for(Int_t adr=0;adr<=47;adr++){
71           Float_t mean  = 150+200*gRandom->Rndm();
72           Float_t sigma = 1+0.3*gRandom->Gaus();
73           if(ddl==7 && row == 10 && dil == 6 ) {mean =  AliHMPIDParam::kPadMeanMasked ; sigma = AliHMPIDParam::kPadSigmaMasked;}                 //to test the dead channel map
74           if(ddl==7 && row == 10 && dil == 4 ) {mean =  AliHMPIDParam::kPadMeanMasked ; sigma = AliHMPIDParam::kPadSigmaMasked;}                 //to test the dead channel map
75           Int_t inhard=((Int_t(mean+nSigmas*sigma))<<9)+Int_t(mean);                                                                             //right calculation, xchecked with Paolo 8/4/2008
76           out << Form("%2i %2i %2i %5.3f %5.3f %4.4x \n",row,dil,adr,mean,sigma,inhard);
77         }
78
79     out.close();
80   }
81   Printf("HMPID - All %i DDL pedestal files created successfully",iDDLmax-iDDLmin+1);
82 }
83 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
84 void SimMap(TMap *pDcsMap,Int_t runTime=1500)
85 {
86   Int_t stepTime=100; //time interval between measurements
87   Int_t startTime=0;
88
89  // phototube current (microAmper)
90  
91   Float_t aArgonCell[] = {2.392914534,1.052269816,0.426667392,0.266823828,0.230182067,0.237871274,0.273670882,0.321822613,0.386860102,0.467150569,0.542086422,0.587703943,0.620307803,0.64317745,0.654298306,0.65143925,0.668053627,0.687757671,0.705103695,0.734520793,0.757498145,0.795069814,0.830174983,0.863526344,0.907381952,0.957449734,0.983424425,1.016224384,1.041821599,1.059507251};
92
93   Float_t aArgonRef[] = {1.149247766,0.496981889,0.186563596,0.10783793,0.085225783,0.084287018,0.095322073,0.109832592,0.130130395,0.157257095,0.180872142,0.196021155,0.205875158,0.212933287,0.217398748,0.218543261,0.223784506,0.233313993,0.243647426,0.254273385,0.270319641,0.288991749,0.310744852,0.333126128,0.359347552,0.383723944,0.412809372,0.432352483,0.455913931,0.477137864}; 
94
95
96   Float_t aFreonCell[] = {0.000104553,0.00011161,0.0002890787,0.0002893622,0.00025613,0.020953063,0.068830326,0.127829805,0.193517566,0.273626387,0.35040611,0.41024375,0.45290783,0.492588729,0.522967756,0.552341878,0.582529366,0.615325809,0.647942483,0.683176041,0.721948683,0.7626279,0.804566145,0.847852349,0.891971886,0.937028468,0.977787018,1.017553926,1.050016046,1.078630805}; 
97
98
99   Float_t aFreonRef[] = {1.125558376,0.475898296,0.178687423,0.100295901,0.078788362,0.077350542,0.087279864,0.101818182,0.121531218,0.147556156,0.170467392,0.184769839,0.194379449,0.202385217,0.209258303,0.215506181,0.222740307,0.231011033,0.241791919,0.253416091,0.2687684,0.286887407,0.310744882,0.332763016,0.356916934,0.383346796,0.409290999,0.433670729,0.454144865,0.474498451};
100
101   TObjArray *pWaveLenght[30], *pArgonCell[30], *pArgonRef[30], *pFreonCell[30], *pFreonRef[30];
102
103   for(Int_t i=0; i<30; i++){
104       
105       pWaveLenght[i] = new TObjArray; pWaveLenght[i]->SetOwner(1);
106       pArgonCell[i]  = new TObjArray; pArgonCell[i]->SetOwner(1);
107       pArgonRef[i]   = new TObjArray; pArgonRef[i]->SetOwner(1);
108       pFreonCell[i]  = new TObjArray; pFreonCell[i]->SetOwner(1);
109       pFreonRef[i]   = new TObjArray; pFreonRef[i]->SetOwner(1);
110
111       pWaveLenght[i]->Add(new AliDCSValue((Float_t)(160+2*i),0));        // wavelenght (nm)
112       pArgonCell[i]->Add(new AliDCSValue((Float_t)(aArgonCell[i]),0));
113       pArgonRef[i]->Add(new AliDCSValue((Float_t)(aArgonRef[i]),0));
114       pFreonRef[i]->Add(new AliDCSValue((Float_t)(aFreonRef[i]),0));
115       pFreonCell[i]->Add(new AliDCSValue((Float_t)(aFreonCell[i]),0));     
116
117
118       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.waveLenght",i)),pWaveLenght[i]);
119       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonCell",i)),pArgonCell[i]);
120       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonReference",i)),pArgonRef[i]);
121       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Cell",i)),pFreonCell[i]);
122       pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Reference",i)),pFreonRef[i]);    
123      }
124   
125   TObjArray *pHV[7];
126   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
127     TObjArray *pPCH4=new TObjArray;  pPCH4->SetOwner(1);
128     TObjArray *pPenv=new TObjArray;  pPenv->SetOwner(1);
129     pHV[iCh]=new TObjArray; pHV[iCh]->SetOwner(1);
130     for(Int_t time=0;time<runTime;time+=stepTime) {
131        pPCH4->Add(new AliDCSValue((Float_t)4.0,time));               //sample CH4 pressure [mBar] differential respect to atm
132        pPenv->Add(new AliDCSValue((Float_t)1000.0 ,time));               //also atm. pressure set to the same value
133 //       pHV[iCh]->Add(new AliDCSValue((Float_t)(1930+iCh*20),time));   //sample chamber HV [V]
134        pHV[iCh]->Add(new AliDCSValue((Float_t)(2050),time));   //sample chamber HV [V]
135     }
136     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
137     pDcsMap->Add(new TObjString(Form("HMP_DET/HMP_ENV/HMP_ENV_PENV.actual.value"                                               )),pPenv);         //atm pressure
138     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
139     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
140     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 
141     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 
142     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 
143     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 
144
145     for(Int_t iRad=0;iRad<3;iRad++){//radiators loop
146       TObjArray *pT1=new TObjArray; pT1->SetOwner(1);
147       TObjArray *pT2=new TObjArray; pT2->SetOwner(1);
148       for (Int_t time=0;time<runTime;time+=stepTime)  pT1->Add(new AliDCSValue((Float_t)(13.0+gRandom->Rndm()),time));  //sample inlet temperature                    Nmean=1.292 @ 13 degrees
149       for (Int_t time=0;time<runTime;time+=stepTime)  pT2->Add(new AliDCSValue((Float_t)(18.0+gRandom->Rndm()),time));  //sample outlet temperature
150       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 
151       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
152     }//radiators loop
153   }//chambers loop
154 }//SimMap()
155 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
156 void DrawInput(TMap *pDcsMap)
157 {
158   
159   TCanvas *c=new TCanvas("cc","Input data",800,800); c->Divide(3,3);
160   
161   AliDCSValue *pVal; Int_t cnt;
162   
163   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
164     if(iCh==6) c->cd(1);  if(iCh==5) c->cd(2);                          
165     if(iCh==4) c->cd(4);  if(iCh==3) c->cd(5);  if(iCh==2) c->cd(6);
166                           if(iCh==1) c->cd(8);  if(iCh==0) c->cd(9); 
167                           
168     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
169     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
170     TGraph *pGrHV=new TGraph; pGrHV->SetMarkerStyle(5); TIter nextHV(pHV);
171     TGraph *pGrP =new TGraph; pGrP ->SetMarkerStyle(5); TIter nextP (pP );
172     
173     for(Int_t iRad=0;iRad<3;iRad++){
174       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);
175       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);
176       TGraph *pGrT1=new TGraph; pGrT1->SetMarkerStyle(5); 
177       TGraph *pGrT2=new TGraph; pGrT2->SetMarkerStyle(5); 
178       pGrT1->SetMaximum(40);
179       pGrT2->SetMaximum(40);
180       cnt=0; while((pVal=(AliDCSValue*)nextT1())) pGrT1->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());
181       cnt=0; while((pVal=(AliDCSValue*)nextT2())) pGrT2->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());
182       if(iRad==0) pGrT1->Draw("AP"); else pGrT1->Draw("same");
183       pGrT2->Draw("same");
184     }//radiators loop
185   }//chambers loop  
186 }//DrawInput()
187 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
188 void DrawOutput()
189 {
190   AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB"); AliCDBManager::Instance()->SetRun(0);
191   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre");
192   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean");
193   AliCDBEntry *pDaqSigEnt=AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");
194   AliCDBEntry *pDaqMaskedEnt=AliCDBManager::Instance()->Get("HMPID/Calib/Masked");
195   AliCDBEntry *pNoiseEnt=AliCDBManager::Instance()->Get("HMPID/Calib/NoiseMap");
196   
197   if(!pQthreEnt || !pNmeanEnt || !pDaqSigEnt || !pDaqMaskedEnt) return;
198   
199   TObjArray *pNmean =(TObjArray*)pNmeanEnt ->GetObject(); 
200   TObjArray *pQthre =(TObjArray*)pQthreEnt ->GetObject(); 
201   TObjArray *pDaqSig=(TObjArray*)pDaqSigEnt->GetObject();
202   TObjArray *pDaqMasked=(TObjArray*)pDaqMaskedEnt->GetObject();
203   TH2F      *pNoiseMap = (TH2F*)pNoiseEnt->GetObject();      
204   
205   TF1 *pRad0,*pRad1,*pRad2;  
206   TCanvas *c2=new TCanvas("c2","Output",800,800); c2->Divide(3,3);
207   TCanvas *c3=new TCanvas("c3","Output Masked",800,800); c3->Divide(3,3);
208   TCanvas *c4=new TCanvas("c4","Noise",800,800); 
209   
210   TH1F *pSig[7];
211   TH2F *phMask[7]; 
212   
213   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
214     TMatrix *pM=(TMatrix*)pDaqSig->At(iCh);
215     TMatrix *pmMask=(TMatrix*)pDaqMasked->At(iCh);
216     
217     pSig[iCh]=new TH1F(Form("sig%i",iCh),"Sigma;[QDC]",100,-5,20); //pSig[iCh]->SetLineColor(iCh+kRed);
218     phMask[iCh]=new TH2F(Form("phMask%i",iCh),"Maksed;pad x;pad y;Q (ADC)",160,0,160,144,0,144); 
219     for(Int_t padx=0;padx<160;padx++) 
220     {
221       for(Int_t pady=0;pady<144;pady++)
222       { 
223         pSig[iCh]->Fill((*pM)(padx,pady));
224         phMask[iCh]->Fill(padx,pady,(*pmMask)(padx,pady));
225        }
226      }
227     c2->cd(7);    if(iCh==0) pSig[iCh]->Draw(); else pSig[iCh]->Draw("same");
228     
229     if(iCh==6) c2->cd(1);  if(iCh==5) c2->cd(2);                          
230     if(iCh==4) c2->cd(4);  if(iCh==3) c2->cd(5);  if(iCh==2) c2->cd(6);
231                            if(iCh==1) c2->cd(8);  if(iCh==0) c2->cd(9); 
232                           
233     TF1 *pRad0=(TF1*)pNmean->At(iCh*3+0); pRad0->Draw();  pRad0->GetXaxis()->SetTimeDisplay(kTRUE); pRad0->GetYaxis()->SetRangeUser(1.28,1.3);
234     TF1 *pRad1=(TF1*)pNmean->At(iCh*3+1); pRad1->Draw("same");
235     TF1 *pRad2=(TF1*)pNmean->At(iCh*3+2); pRad2->Draw("same");
236      
237     if(iCh==6) c3->cd(1);  if(iCh==5) c3->cd(2);                          
238     if(iCh==4) c3->cd(4);  if(iCh==3) c3->cd(5);  if(iCh==2) c3->cd(6);
239                            if(iCh==1) c3->cd(8);  if(iCh==0) c3->cd(9); 
240      phMask[iCh]->Draw("colz");
241     
242   }//chambers loop   
243   
244   c4->cd();
245   pNoiseMap->Draw("colz");
246   
247    
248 }//DrawOutput()
249 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++