1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 /*-------------------------------------------------------------------------------
19 /* 14/12/07 New version: MUONTRKda.cxx,v 1.10
20 /*-------------------------------------------------------------------------------
22 Version for MUONTRKda MUON tracking
23 (A. Baldisseri, J.-L. Charvet & Ch. Finck)
26 Rem: AliMUON2DMap stores all channels, even those which are not connected
27 if pedMean == -1, channel not connected to a pad
38 #include <Riostream.h>
43 #include "AliMUONLogger.h"
44 #include "AliMUONRawStreamTracker.h"
45 #include "AliMUONDspHeader.h"
46 #include "AliMUONBlockHeader.h"
47 #include "AliMUONBusStruct.h"
48 #include "AliMUONDDLTracker.h"
49 #include "AliMUONVStore.h"
50 #include "AliMUON2DMap.h"
51 #include "AliMUONCalibParamND.h"
52 #include "AliMpIntPair.h"
53 #include "AliMpConstants.h"
54 #include "AliRawReaderDate.h"
62 #include "TStopwatch.h"
64 #include "TTimeStamp.h"
65 #include "TGraphErrors.h"
68 #include "TPluginManager.h"
74 const Int_t gkNChannels = AliMpConstants::ManuNofChannels();
75 const Int_t gkADCMax = 4095;
77 AliMUONVStore* gPedestalStore = new AliMUON2DMap(kFALSE);
81 UInt_t gRunNumber = 0;
83 Int_t gNDateEvents = 0;
84 Int_t gPrintLevel = 1; // global printout variable (others: 2 and 3)
85 Int_t gPlotLevel = 0; // global plot variable
87 TH1F* gPedMeanHisto = 0x0;
88 TH1F* gPedSigmaHisto = 0x0;
89 Char_t gHistoFileName[256];
92 Char_t gHistoFileName_gain[256]="MUONTRKda_gain_data.root";
93 Char_t gRootFileName[256];
94 Char_t gOutFolder[256]=".";
96 Char_t filenam[256]="MUONTRKda_gain";
101 TString flatOutputFile;
102 TString logOutputFile;
103 TString gCommand("ped");
110 Double_t funcLin(Double_t *x, Double_t *par)
112 return par[0] + par[1]*x[0];
116 Double_t funcParabolic(Double_t *x, Double_t *par)
118 return par[0]*x[0]*x[0];
122 Double_t funcCalib(Double_t *x, Double_t *par)
124 Double_t xLim= par[3];
126 if(x[0] <= xLim) return par[0] + par[1]*x[0];
128 Double_t yLim = par[0]+ par[1]*xLim;
129 return yLim + par[1]*(x[0] - xLim) + par[2]*(x[0] - xLim)*(x[0] - xLim);
134 void MakePed(Int_t busPatchId, Int_t manuId, Int_t channelId, Int_t charge)
137 AliMUONVCalibParam* ped =
138 static_cast<AliMUONVCalibParam*>(gPedestalStore->FindObject(busPatchId, manuId));
142 ped = new AliMUONCalibParamND(2, gkNChannels,busPatchId, manuId, -1.); // put default wise -1, not connected channel
143 gPedestalStore->Add(ped);
147 ped->SetValueAsDouble(channelId, 0, 0.);
148 ped->SetValueAsDouble(channelId, 1, 0.);
151 Double_t pedMean = ped->ValueAsDouble(channelId, 0) + charge;
152 Double_t pedSigma = ped->ValueAsDouble(channelId, 1) + charge*charge;
154 ped->SetValueAsDouble(channelId, 0, pedMean);
155 ped->SetValueAsDouble(channelId, 1, pedSigma);
160 void MakePedStore(TString flatOutputFile = "")
171 TFile* histoFile = 0;
173 if (gCommand.CompareTo("ped") == 0)
175 sprintf(gHistoFileName,"%s/MUONTRKda_ped_%d.root",gOutFolder,gRunNumber);
176 histoFile = new TFile(gHistoFileName,"RECREATE","MUON Tracking pedestals");
180 sprintf(name,"pedmean_allch");
181 sprintf(title,"Pedestal mean all channels");
185 gPedMeanHisto = new TH1F(name,title,nx,xmin,xmax);
186 gPedMeanHisto->SetDirectory(histoFile);
188 sprintf(name,"pedsigma_allch");
189 sprintf(title,"Pedestal sigma all channels");
193 gPedSigmaHisto = new TH1F(name,title,nx,xmin,xmax);
194 gPedSigmaHisto->SetDirectory(histoFile);
196 tree = new TTree("t","Pedestal tree");
197 tree->Branch("bp",&busPatchId,"bp/I");
198 tree->Branch("manu",&manuId,",manu/I");
199 tree->Branch("channel",&channelId,",channel/I");
200 tree->Branch("pedMean",&pedMean,",pedMean/D");
201 tree->Branch("pedSigma",&pedSigma,",pedSigma/D");
204 if (!flatOutputFile.IsNull()) {
205 fileout.open(flatOutputFile.Data());
206 fileout<<"//===========================================================================" << endl;
207 fileout<<"// Pedestal file calculated by MUONTRKda"<<endl;
208 fileout<<"//===========================================================================" << endl;
209 fileout<<"// * Run : " << gRunNumber << endl;
210 fileout<<"// * Date : " << date.AsString("l") <<endl;
211 fileout<<"// * Statictics : " << gNEvents << endl;
212 fileout<<"// * # of MANUS : " << gNManu << endl;
213 fileout<<"// * # of channels : " << gNChannel << endl;
215 fileout<<"//---------------------------------------------------------------------------" << endl;
216 fileout<<"//---------------------------------------------------------------------------" << endl;
217 fileout<<"// BP MANU CH. MEAN SIGMA"<<endl;
218 fileout<<"//---------------------------------------------------------------------------" << endl;
222 // iterator over pedestal
223 TIter next(gPedestalStore->CreateIterator());
224 AliMUONVCalibParam* ped;
226 while ( ( ped = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
228 busPatchId = ped->ID0();
231 for (channelId = 0; channelId < ped->Size() ; ++channelId) {
233 pedMean = ped->ValueAsDouble(channelId, 0);
235 if (pedMean > 0) { // connected channels
237 ped->SetValueAsDouble(channelId, 0, pedMean/(Double_t)gNEvents);
239 pedMean = ped->ValueAsDouble(channelId, 0);
240 pedSigma = ped->ValueAsDouble(channelId, 1);
242 ped->SetValueAsDouble(channelId, 1, TMath::Sqrt(TMath::Abs(pedSigma/(Double_t)gNEvents - pedMean*pedMean)));
244 pedMean = ped->ValueAsDouble(channelId, 0) + 0.5 ;
245 // pedMean = ped->ValueAsDouble(channelId, 0) ;
246 pedSigma = ped->ValueAsDouble(channelId, 1);
249 if (!flatOutputFile.IsNull()) {
250 fileout << "\t" << busPatchId << "\t" << manuId <<"\t"<< channelId << "\t"
251 << pedMean <<"\t"<< pedSigma << endl;
254 if (gCommand.CompareTo("ped") == 0)
256 gPedMeanHisto->Fill(pedMean);
257 gPedSigmaHisto->Fill(pedSigma);
266 if (!flatOutputFile.IsNull()) fileout.close();
268 if (gCommand.CompareTo("ped") == 0)
279 void MakePedStoreForGain(Int_t injCharge)
281 // store pedestal map in root file
283 // Int_t injCharge = 200;
287 // compute and store pedestals
288 sprintf(flatFile,"%s/%s_%d_DAC_%d.ped",gOutFolder,filenam,gRunNumber,injCharge);
289 cout << "\nMUONTRKda : Flat file generated : " << flatFile << "\n";
290 MakePedStore(flatFile);
292 TString mode("UPDATE");
294 if (gCommand.Contains("cre")) {
297 TFile* histoFile = new TFile(gHistoFileName_gain, mode.Data(), "MUON Tracking Gains");
299 // second argument should be the injected charge, taken from config crocus file
300 // put also info about run number could be usefull
301 AliMpIntPair* pair = new AliMpIntPair(gRunNumber, injCharge);
303 if (mode.CompareTo("UPDATE") == 0) {
304 tree = (TTree*)histoFile->Get("t");
305 tree->SetBranchAddress("run",&pair);
306 tree->SetBranchAddress("ped",&gPedestalStore);
309 tree = new TTree("t","Pedestal tree");
310 tree->Branch("run", "AliMpIntPair",&pair);
311 tree->Branch("ped", "AliMUON2DMap",&gPedestalStore);
312 tree->SetBranchAddress("run",&pair);
313 tree->SetBranchAddress("ped",&gPedestalStore);
318 tree->Write("t", TObject::kOverwrite); // overwrite the tree
325 // void MakeGainStore(TString flatOutputFile)
328 Double_t goodA1Min = 0.5;
329 Double_t goodA1Max = 2.;
330 Double_t goodA2Min = -0.5E-03;
331 Double_t goodA2Max = 1.E-03;
333 // open file mutrkgain.root
334 // read again the pedestal for the calibration runs (9 runs ?)
335 // need the injection charge from config file (to be done)
336 // For each channel make a TGraphErrors (mean, sigma) vs injected charge
337 // Fit with a polynomial fct
338 // store the result in a flat file.
341 TFile* histoFile = new TFile(gHistoFileName_gain);
343 AliMUON2DMap* map[11];
344 AliMUONVCalibParam* ped[11];
345 AliMpIntPair* run[11];
347 //read back from root file
348 TTree* tree = (TTree*)histoFile->Get("t");
349 Int_t nEntries = tree->GetEntries();
352 for (Int_t i = 0; i < nEntries; ++i) {
355 tree->SetBranchAddress("ped",&map[i]);
356 tree->SetBranchAddress("run",&run[i]);
358 // std::cout << map[i] << " " << run[i] << std::endl;
360 gRunNumber=(UInt_t)run[0]->GetFirst();
363 cout<<"\n ******** MUONTRKda for Gain computing (Run = " << gRunNumber << ")\n" << endl;
364 cout<<" * Date : " << date.AsString("l") << "\n" << endl;
365 cout << " Entries = " << nEntries << " DAC values \n" << endl;
366 for (Int_t i = 0; i < nEntries; ++i) {
367 cout<< " Run = " << (Double_t)run[i]->GetFirst() << " DAC = " << (Double_t)run[i]->GetSecond() << endl;
372 Double_t pedMean[11];
373 Double_t pedSigma[11];
374 Double_t injCharge[11];
375 Double_t injChargeErr[11];
379 sprintf(filename,"%s/%s_%d.log",gOutFolder,filenam,gRunNumber);
380 logOutputFile=filename;
382 filcout.open(logOutputFile.Data());
383 filcout<<"//====================================================" << endl;
384 filcout<<"// MUONTRKda for Gain computing (Run = " << gRunNumber << ")" << endl;
385 filcout<<"//====================================================" << endl;
386 filcout<<"// * Date : " << date.AsString("l") << "\n" << endl;
391 // why 2 files ? (Ch. F.)
396 sprintf(filename,"%s/%s_%d.param",gOutFolder,filenam,gRunNumber);
397 cout << " fit parameter file = " << filename << "\n";
398 pfilen = fopen (filename,"w");
400 fprintf(pfilen,"//===================================================================\n");
401 fprintf(pfilen,"// BP MANU CH. a0 a1 a2 xlim P(chi2) P(chi2)2 Q\n");
402 fprintf(pfilen,"//===================================================================\n");
403 fprintf(pfilen,"// * Run : %d \n",gRunNumber);
404 fprintf(pfilen,"//===================================================================\n");
406 sprintf(filename,"%s/%s_%d.bad",gOutFolder,filenam,gRunNumber);
407 cout << " Bad channel file = " << filename << "\n";
408 pfilef = fopen (filename,"w");
410 fprintf(pfilef,"//=================================================\n");
411 fprintf(pfilef,"// Bad Channel file calculated by MUONTRKda \n");
412 fprintf(pfilef,"//=================================================\n");
413 fprintf(pfilef,"// * Run : %d \n",gRunNumber);
414 fprintf(pfilef,"// * Date : %s \n",date.AsString("l"));
415 fprintf(pfilef,"//=======================================\n");
416 fprintf(pfilef,"// BP MANU CH. a1 a2 thres. Q\n");
417 fprintf(pfilef,"//=======================================\n");
421 if(flatOutputFile.IsNull())
423 sprintf(filename,"%s_%d.par",filenam,gRunNumber);
424 flatOutputFile=filename;
426 if(!flatOutputFile.IsNull())
428 pfilew = fopen (flatOutputFile.Data(),"w");
430 fprintf(pfilew,"//=================================================\n");
431 fprintf(pfilew,"// Calibration file calculated by MUONTRKda \n");
432 fprintf(pfilew,"//=================================================\n");
433 fprintf(pfilew,"// * Run : %d \n",gRunNumber);
434 fprintf(pfilew,"// * Date : %s \n",date.AsString("l"));
435 fprintf(pfilew,"// * Statictics : %d \n",gNEvents);
436 fprintf(pfilew,"// * # of MANUS : %d \n",gNManu);
437 fprintf(pfilew,"// * # of channels : %d \n",gNChannel);
438 fprintf(pfilew,"//-------------------------------------------------\n");
439 fprintf(pfilew,"//=======================================\n");
440 fprintf(pfilew,"// BP MANU CH. a1 a2 thres. Q\n");
441 fprintf(pfilew,"//=======================================\n");
447 sprintf(filename,"%s/%s_%d.peak",gOutFolder,filenam,gRunNumber);
448 cout << " File containing Peak mean values = " << filename << "\n";
449 pfilep = fopen (filename,"w");
451 fprintf(pfilep,"//===============================================================================================================================\n");
452 fprintf(pfilep,"// * Run : %d \n",gRunNumber);
453 fprintf(pfilep,"//===============================================================================================================================\n");
454 fprintf(pfilep,"// BP MANU CH. Ped. <0> <1> <2> <3> <4> <5> <6> <7> <8> <9> <10> \n");
455 // fprintf(pfilep,"// %9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f fC\n",level_fC[0],level_fC[1],level_fC[2],level_fC[3],level_fC[4],level_fC[5],level_fC[6],level_fC[7],level_fC[8],level_fC[9],level_fC[10]);
456 // fprintf(pfilep,"// %9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f\n",level_err[0],level_err[1],level_err[2],level_err[3],level_err[4],level_err[5],level_err[6],level_err[7],level_err[8],level_err[9],level_err[10]);
457 fprintf(pfilep,"//===============================================================================================================================\n");
464 TFile* gainFile = 0x0;
465 sprintf(gRootFileName,"%s/%s_%d.root",gOutFolder,filenam,gRunNumber);
466 gainFile = new TFile(gRootFileName,"RECREATE");
469 Double_t chi2P2 = 0.;
471 Double_t prChi2P2 =0;
481 Double_t capa=0.2; // internal capacitor (pF)
483 TTree *tg = new TTree("tg","TTree avec class Manu_DiMu");
485 tg->Branch("bp",&busPatchId, "busPatchId/I");
486 tg->Branch("manu",&manuId, "manuId/I");
487 tg->Branch("channel",&channelId, "channelId/I");
489 tg->Branch("a0",&a0, "a0/D");
490 tg->Branch("a1",&a1, "a1/D");
491 tg->Branch("a2",&a2, "a2/D");
492 tg->Branch("Pchi2",&prChi2, "prChi2/D");
493 tg->Branch("Pchi2_2",&prChi2P2, "prChi2P2/D");
494 tg->Branch("Threshold",&threshold, "threshold/I");
495 tg->Branch("Q",&Q, "Q/I");
496 tg->Branch("p1",&p1, "p1/I");
497 tg->Branch("p2",&p2, "p2/I");
498 tg->Branch("gain",&gain, "gain/D");
500 // bad BusPatch and manu
501 Int_t num_tot_BP=800;
502 Int_t num_tot_Manu=1500;
503 // Int_t bad_channel[num_tot_BP][num_tot_Manu];
504 Int_t bad_channel[800][1500];
505 for ( Int_t i = 0; i < num_tot_BP ; i++ )
506 { for ( Int_t j = 0; j < num_tot_Manu ; j++ ) bad_channel[i][j]=0;}
511 // iterates over the first pedestal run
512 TIter next(map[0]->CreateIterator());
513 AliMUONVCalibParam* p;
516 Int_t nGoodChannel = 0;
517 Int_t nGoodChannel_a1 = 0;
518 Int_t nBadChannel = 0;
519 Int_t nBadChannel_a1 = 0;
520 Int_t nBadChannel_a2 = 0;
522 Double_t sumProbChi2 = 0.;
524 Double_t sumProbChi2P2 = 0.;
527 Double_t x[11], xErr[11], y[11], yErr[11];
528 Double_t xp[11], xpErr[11], yp[11], ypErr[11];
530 while ( ( p = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
534 busPatchId = p->ID0();
537 // read back pedestal from the other runs for the given (bupatch, manu)
538 for (Int_t i = 1; i < nEntries; ++i) {
539 ped[i] = static_cast<AliMUONVCalibParam*>(map[i]->FindObject(busPatchId, manuId));
542 // compute for each channel the gain parameters
543 for ( channelId = 0; channelId < ped[0]->Size() ; ++channelId ) {
548 for (Int_t i = 0; i < nEntries; ++i) {
550 if (!ped[i]) continue; //shouldn't happen.
551 pedMean[i] = ped[i]->ValueAsDouble(channelId, 0);
552 pedSigma[i] = ped[i]->ValueAsDouble(channelId, 1);
553 injCharge[i] = (Double_t)run[i]->GetSecond();
554 injChargeErr[i] = 0.01*injCharge[i];
555 if(injChargeErr[i] <= 1.) injChargeErr[i]=1.;
557 if (pedMean[i] < 0) continue; // not connected
559 if (pedSigma[i] <= 0) pedSigma[i] = 1.; // should not happen.
564 // print_peak_mean_values
568 fprintf(pfilep,"%4i%5i%5i%10.3f",busPatchId,manuId,channelId,0.);
569 fprintf(pfilep,"%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f mV\n",pedMean[0],pedMean[1],pedMean[2],pedMean[3],pedMean[4],pedMean[5],pedMean[6],pedMean[7],pedMean[8],pedMean[9],pedMean[10]);
570 fprintf(pfilep," %9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f%9.1f \n",pedSigma[0],pedSigma[1],pedSigma[2],pedSigma[3],pedSigma[4],pedSigma[5],pedSigma[6],pedSigma[7],pedSigma[8],pedSigma[9],pedSigma[10]);
576 // Fit Method: Linear fit over 6 points + fit parabolic function over 3 points)
578 // 1. - linear fit over 6 points
580 Double_t par[4] = {0.,0.,0.,gkADCMax};
583 Int_t nbs = nEntries - nInit;
584 Int_t nbpf1 = 6; // linear fit over nbf1 points
586 for (Int_t j = 0; j < nbs; ++j)
590 xErr[j] = pedSigma[k];
592 yErr[j] = injChargeErr[k];
596 TF1 *f1 = new TF1("f1",funcLin,0.,gkADCMax,2);
597 TGraphErrors *graphErr = new TGraphErrors(nbpf1, x, y, xErr, yErr);
599 f1->SetParameters(0,0);
601 graphErr->Fit("f1","RQ");
603 chi2 = f1->GetChisquare();
604 f1->GetParameters(par);
610 prChi2 = TMath::Prob(chi2, nbpf1 - 2);
612 Double_t xLim = pedMean[nInit + nbpf1 - 1];
613 Double_t yLim = par[0]+par[1] * xLim;
618 // 2. - Translation : new origin (xLim, yLim) + parabolic fit over nbf2 points
620 Int_t nbpf2 = nEntries - (nInit + nbpf1) + 1;
624 for (Int_t j = 0; j < nbpf2; j++)
626 Int_t k = j + (nInit + nbpf1) - 1;
627 xp[j] = pedMean[k] - xLim;
628 xpErr[j] = pedSigma[k];
630 yp[j] = injCharge[k] - yLim - par[1]*xp[j];
631 ypErr[j] = injChargeErr[k];
635 TF1 *f2 = new TF1("f2",funcParabolic,0.,gkADCMax,1);
636 TGraphErrors *graphErr = new TGraphErrors(nbpf2, xp, yp, xpErr, ypErr);
638 graphErr->Fit(f2,"RQ");
639 chi2P2 = f2->GetChisquare();
640 f2->GetParameters(par);
646 prChi2P2 = TMath::Prob(chi2P2, nbpf2-1);
649 // ------------- print out in log file
650 // if (busPatchId == 6 && manuId == 116 && ( channelId >= 17 && channelId <= 20) )
652 // filcout << " \n ********! Print_out.: BP= " << busPatchId << " Manu_Id= " << manuId
653 // << " Ch.= " << channelId << ":" << endl;
655 // for (Int_t j = 0; j < nbpf1; ++j)
656 // {filcout << j << " " << x[j] << " " << xErr[j] << " " << y[j] << " " << yErr[j] << endl;}
657 // filcout << " a0,a1 = " << a0 << " , " << a1 << " pr_chi2 = " << prChi2 << endl ;
659 // for (Int_t j = 0; j < nbpf2; ++j)
660 // {filcout << j << " " << xp[j] << " " << xpErr[j] << " " << yp[j] << " " << ypErr[j] << endl;}
661 // filcout << " a2 = " << par[0] << " pr_chi2_2 = " << prChi2P2 << endl;
664 // ------------------------------------------
681 p1 = TMath::Nint(ceil(prChi2*14))+1;
682 p2 = TMath::Nint(ceil(prChi2P2*14))+1;
684 Double_t x0 = -par[0]/par[1]; // value of x corresponding to à 0 fC
685 threshold = TMath::Nint(ceil(par[3]-x0)); // linear if x < threshold
689 fprintf(pfilen,"%4i %4i %2i",busPatchId,manuId,channelId);
690 fprintf(pfilen," %6.2f %6.4f %10.3e %4.2f %5.3f %x %5.3f %x\n",
691 par[0], par[1], par[2], par[3], prChi2, p1, prChi2P2, p2);
696 if(par[1]< goodA1Min || par[1]> goodA1Max)
700 if (gPrintLevel && nBadChannel_a1 < 1)
702 cout << " !!!!! " << nBadChannel_a1 << " !!!!!!!! Bad Calib.: BP= " << busPatchId << " Manu_Id= " << manuId <<
703 " Ch.= " << channelId << ":";
704 cout << " a1 = " << par[1] << " out of limit : [" << goodA1Min << "," << goodA1Max <<
709 if(par[2]< goodA2Min || par[2]> goodA2Max)
713 if (gPrintLevel && nBadChannel_a2 < 1)
715 cout << " !!!!! " << nBadChannel_a2 << " !!!!!!!! Bad Calib.: BP= " << busPatchId << " Manu_Id= " << manuId
716 << " Ch.= " << channelId << ":";
717 cout << " a2 = " << par[2] << " out of limit : [" << goodA2Min << "," << goodA2Max
720 for (Int_t j = 0; j < nbpf2; ++j)
721 {cout << j << " " << x[j] << " " << xErr[j] << " " << y[j] << " " << yErr[j] << endl;}
725 Q = p1*16 + p2; // fit quality
726 if(p1==0)Q=0; // bad linear fit <=> bad calibration
731 sumProbChi2P2 += prChi2P2;
737 if(busPatchId < num_tot_BP && manuId < num_tot_Manu) bad_channel[busPatchId][manuId]++;
738 else{cout << " Warning : busPatch = " << busPatchId << " Manu = " << manuId << endl;}
739 if(gPrintLevel>=2)fprintf(pfilef,"%4i %5i %2i %7.4f %10.3e %4i %2x\n",busPatchId,manuId,channelId,par[1],par[2],threshold,Q);
746 sumProbChi2 += prChi2;
748 gain=1./(par[1]*capa);
754 if (!flatOutputFile.IsNull())
756 fprintf(pfilew,"%4i %5i %2i %7.4f %10.3e %4i %2x\n",busPatchId,manuId,channelId,par[1],par[2],threshold,Q);
762 if(Q==0 and nplot < 100)
763 // if(p1>1 && p2==0 and nplot < 100)
764 // if(p1>1 && p2>1 and nplot < 100)
767 TF1 *f2Calib = new TF1("f2Calib",funcCalib,0.,gkADCMax,NFITPARAMS);
769 TGraphErrors *graphErr = new TGraphErrors(nEntries,pedMean,injCharge,pedSigma,injChargeErr);
771 sprintf(graphName,"BusPatch_%d_Manu_%d_Ch_%d",busPatchId, manuId,channelId);
773 graphErr->SetTitle(graphName);
774 graphErr->SetMarkerColor(3);
775 graphErr->SetMarkerStyle(12);
776 graphErr->Write(graphName);
778 sprintf(graphName,"f2_BusPatch_%d_Manu_%d_Ch_%d",busPatchId, manuId,channelId);
779 f2Calib->SetTitle(graphName);
780 f2Calib->SetLineColor(4);
781 f2Calib->SetParameters(par);
782 f2Calib->Write(graphName);
790 if(fmod(nmanu,100)==0)std::cout << " Nb manu = " << nmanu << std::endl;
793 // file outputs for gain
794 if (!flatOutputFile.IsNull()) fclose(pfilew);
795 if(gPrintLevel==2){ fclose(pfilen); fclose(pfilef);}
796 if(gPrintLevel==3) fclose(pfilep);
804 filcout << "\n List of problematic BusPatch and Manu " << endl;
805 filcout << " ========================================" << endl;
806 filcout << " BP Manu Nb Channel " << endl ;
807 filcout << " ========================================" << endl;
808 for ( Int_t i = 0 ; i < num_tot_BP ; i++ )
809 { for ( Int_t j = 0 ; j < num_tot_Manu ; j++ )
810 if (bad_channel[i][j] != 0 ) filcout << "\t" << i << "\t " << j << "\t\t" << bad_channel[i][j] << endl;}
811 filcout << " ========================================" << endl;
814 filcout << "\n Nb of channels in raw data = " << nmanu*64 << " (" << nmanu << " Manu)" << endl;
815 filcout << "\n Nb of fully calibrated channel = " << nGoodChannel << " (" << goodA1Min << "<a1<" << goodA1Max
816 << " and " << goodA2Min << "<a2<" << goodA2Max << ") " << endl;
817 filcout << "\n Nb of Bad channel = " << nBadChannel << endl;
819 filcout << "\n Nb of Good a1 channels = " << nGoodChannel_a1 << " (" << goodA1Min << "<a1<" << goodA1Max << ") " << endl;
821 Double_t meanA1 = sumA1/(nGoodChannel_a1);
822 Double_t meanProbChi2 = sumProbChi2/(nGoodChannel_a1);
823 Double_t meanA2 = sumA2/(nGoodChannel);
824 Double_t meanProbChi2P2 = sumProbChi2P2/(nGoodChannel);
826 Double_t capaManu = 0.2; // pF
827 filcout << "\n linear fit : <a1> = " << meanA1 << "\t <gain> = " << 1./(meanA1*capaManu)
828 << " mV/fC (capa= " << capaManu << " pF)" << endl;
829 filcout << " Prob(chi2)> = " << meanProbChi2 << endl;
830 filcout << "\n parabolic fit: <a2> = " << meanA2 << endl;
831 filcout << " Prob(chi2)> = " << meanProbChi2P2 << "\n" << endl;
840 //*************************************************************//
843 int main(Int_t argc, Char_t **argv)
846 // needed for streamer application
847 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
853 TFitter *minuitFit = new TFitter(NFITPARAMS);
854 TVirtualFitter::SetFitter(minuitFit);
856 Int_t skipEvents = 0;
857 Int_t maxEvents = 1000000;
858 Int_t MaxDateEvents = 1000000;
861 Int_t threshold = -1;
862 Char_t inputFile[256];
864 Int_t gGlitchErrors= 0;
865 Int_t gParityErrors= 0;
866 Int_t gPaddingErrors= 0;
868 TString fesOutputFile;
869 TString crocusOutputFile;
870 TString crocusConfigFile;
874 // decode the input line
875 for (Int_t i = 1; i < argc; i++) // argument 0 is the executable name
880 if (arg[0] != '-') continue;
885 sprintf(inputFile,argv[i]);
889 flatOutputFile = argv[i];
893 sprintf(gOutFolder,argv[i]);
897 crocusOutputFile = argv[i];
901 crocusConfigFile = argv[i];
909 gPrintLevel=atoi(argv[i]);
913 gPlotLevel=atoi(argv[i]);
917 skipEvents=atoi(argv[i]);
921 injCharge=atoi(argv[i]);
925 sscanf(argv[i],"%d",&MaxDateEvents);
929 sscanf(argv[i],"%d",&maxEvents);
933 sscanf(argv[i],"%lf",&nSigma);
937 sprintf(gHistoFileName_gain,argv[i]);
941 sscanf(argv[i],"%d",&threshold);
945 printf("\n******************* %s usage **********************",argv[0]);
946 printf("\n%s -options, the available options are :",argv[0]);
947 printf("\n-h help (this screen)");
950 printf("\n-f <raw data file> (default = %s)",inputFile);
951 printf("\n-c <Crocus config. file> (default = %s)",crocusConfigFile.Data());
954 printf("\n-a <Flat ASCII file> (default = %s)",flatOutputFile.Data());
955 printf("\n-o <CROCUS cmd file> (default = %s)",crocusOutputFile.Data());
957 printf("\n Options");
958 printf("\n-b <output directory> (default = %s)",gOutFolder);
959 printf("\n-d <print level> (default = %d)",gPrintLevel);
960 printf("\n-g <plot level> (default = %d)",gPlotLevel);
961 printf("\n-l <DAC level> (default = %d)",injCharge);
962 printf("\n-m <max date events> (default = %d)",MaxDateEvents);
963 printf("\n-s <skip events> (default = %d)",skipEvents);
964 printf("\n-n <max events> (default = %d)",maxEvents);
965 printf("\n-p <n sigmas> (default = %f)",nSigma);
966 printf("\n-r root file data for gain(default = %s)",gHistoFileName_gain);
967 printf("\n-t <threshold (-1 = no)> (default = %d)",threshold);
968 printf("\n-e <execute ped/gain> (default = %s)",gCommand.Data());
969 printf("\n-e <gain create> make gain & create a new root file");
970 printf("\n-e <gain> make gain & update root file");
971 printf("\n-e <gain compute> make gain & compute gains");
976 printf("%s : bad argument %s (please check %s -h)\n",argv[0],argv[i],argv[0]);
977 argc = 2; exit(-1); // exit if error
981 // set gCommand to lower case
985 // decoding the events
991 gPedSigmaHisto = 0x0;
997 // once we have a configuration file in db
998 // copy locally a file from daq detector config db
999 // The current detector is identified by detector code in variable
1000 // DATE_DETECTOR_CODE. It must be defined.
1001 // If environment variable DAQDA_TEST_DIR is defined, files are copied from DAQDA_TEST_DIR
1002 // instead of the database. The usual environment variables are not needed.
1003 if (!crocusConfigFile.IsNull()) {
1004 status = daqDA_DB_getFile("myconfig", crocusConfigFile.Data());
1006 printf("Failed to get config file : %d\n",status);
1012 status = monitorSetDataSource(inputFile);
1014 cerr << "ERROR : monitorSetDataSource status (hex) = " << hex << status
1015 << " " << monitorDecodeError(status) << endl;
1018 status = monitorDeclareMp("MUON Tracking monitoring");
1020 cerr << "ERROR : monitorDeclareMp status (hex) = " << hex << status
1021 << " " << monitorDecodeError(status) << endl;
1029 TString key("MUONTRKda :");
1031 AliMUONRawStreamTracker* rawStream = 0;
1034 if (gCommand.CompareTo("comp") != 0)
1036 cout << "\nMUONTRKda : Reading data from file " << inputFile <<endl;
1040 if (gNDateEvents >= MaxDateEvents) break;
1041 if (gNEvents >= maxEvents) break;
1042 if (gNEvents && gNEvents % 100 == 0)
1043 cout<<"Cumulated: DATE events = " << gNDateEvents << " Used events = " << gNEvents << endl;
1045 // check shutdown condition
1046 if (daqDA_checkShutdown())
1049 // Skip Events if needed
1050 while (skipEvents) {
1051 status = monitorGetEventDynamic(&event);
1056 status = monitorGetEventDynamic(&event);
1058 cout<<"EOF found"<<endl;
1062 // decoding rawdata headers
1063 AliRawReader *rawReader = new AliRawReaderDate(event);
1065 Int_t eventType = rawReader->GetType();
1066 gRunNumber = rawReader->GetRunNumber();
1068 // Output log file initialisations
1072 if (gCommand.CompareTo("ped") == 0){
1073 sprintf(flatFile,"%s/MUONTRKda_ped_%d.log",gOutFolder,gRunNumber);
1074 logOutputFile=flatFile;
1076 filcout.open(logOutputFile.Data());
1077 filcout<<"//=================================================" << endl;
1078 filcout<<"// MUONTRKda for Pedestal run = " << gRunNumber << endl;
1079 cout<<"\n ******** MUONTRKda for Pedestal run = " << gRunNumber << "\n" << endl;
1082 if (gCommand.Contains("gain")){
1083 sprintf(flatFile,"%s/%s_%d_DAC_%d.log",gOutFolder,filenam,gRunNumber,injCharge);
1084 logOutputFile=flatFile;
1086 filcout.open(logOutputFile.Data());
1087 filcout<<"//=================================================" << endl;
1088 filcout<<"// MUONTRKda for Gain run = " << gRunNumber << " (DAC=" << injCharge << ")" << endl;
1089 cout<<"\n ******** MUONTRKda for Gain run = " << gRunNumber << " (DAC=" << injCharge << ")\n" << endl;
1092 filcout<<"//=================================================" << endl;
1093 filcout<<"// * Date : " << date.AsString("l") << "\n" << endl;
1094 cout<<" * Date : " << date.AsString("l") << "\n" << endl;
1102 if (eventType != PHYSICS_EVENT)
1103 continue; // for the moment
1105 // decoding MUON payload
1106 // AliMUONRawStreamTracker* rawStream = new AliMUONRawStreamTracker(rawReader);
1107 rawStream = new AliMUONRawStreamTracker(rawReader);
1108 rawStream->DisableWarnings();
1110 // // loops over DDL
1111 // rawStream->First(); // if GlitchError ? what we are doing ?
1112 // while( (status = rawStream->Next(busPatchId, manuId, channelId, charge)) )
1115 // if (gNEvents == 0) gNChannel++;
1117 // MakePed(busPatchId, (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
1121 // if (!rawStream->IsErrorMessage()) {
1125 // loops over DDL to find good events (Alberto 11/12/07)
1126 rawStream->First(); // if GlitchError ? what we are doing ?
1127 while( (status = rawStream->Next(busPatchId, manuId, channelId, charge)) ) {
1130 if (!rawStream->IsErrorMessage()) {
1131 // loops over DDL to find good events
1132 rawStream->First(); // if GlitchError ? what we are doing ?
1133 while( (status = rawStream->Next(busPatchId, manuId, channelId, charge)) ) {
1138 MakePed(busPatchId, (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
1144 filcout<<"Event # "<<*(rawReader->GetEventId())<<" rejected"<<endl;
1146 if (rawStream->GetPayLoad()->GetGlitchErrors()) gGlitchErrors++;
1147 if (rawStream->GetPayLoad()->GetParityErrors()) gParityErrors++;
1148 if (rawStream->GetPayLoad()->GetPaddingErrors()) gPaddingErrors++;
1150 AliMUONLogger* log = rawStream->GetPayLoad()->GetErrorLogger();
1151 log->Print(key, filcout);
1161 if (gCommand.CompareTo("ped") == 0)
1163 sprintf(flatFile,"MUONTRKda_ped_%d.ped",gRunNumber);
1164 if(flatOutputFile.IsNull())flatOutputFile=flatFile;
1165 MakePedStore(flatOutputFile);
1168 // option gain -> update root file with pedestal results
1169 // gain + create -> recreate root file
1170 // gain + comp -> update root file and compute gain parameters
1172 if (gCommand.Contains("gain"))
1174 MakePedStoreForGain(injCharge);
1177 if (gCommand.Contains("comp"))
1180 // if(flatOutputFile.IsNull())flatOutputFile="MUONTRKda_gain.par";
1181 // MakeGainStore(flatOutputFile);
1186 delete gPedestalStore;
1189 TVirtualFitter::SetFitter(0);
1193 if (gCommand.CompareTo("comp") != 0)
1195 cout << "\nMUONTRKda : Nb of DATE events = " << gNDateEvents << endl;
1196 cout << "MUONTRKda : Nb of Glitch errors = " << gGlitchErrors << endl;
1197 cout << "MUONTRKda : Nb of Parity errors = " << gParityErrors << endl;
1198 cout << "MUONTRKda : Nb of Padding errors = " << gPaddingErrors << endl;
1199 cout << "MUONTRKda : Nb of events used = " << gNEvents << endl;
1201 filcout << "\nMUONTRKda : Nb of DATE events = " << gNDateEvents << endl;
1202 filcout << "MUONTRKda : Nb of Glitch errors = " << gGlitchErrors << endl;
1203 filcout << "MUONTRKda : Nb of Parity errors = " << gParityErrors << endl;
1204 filcout << "MUONTRKda : Nb of Padding errors = " << gPaddingErrors << endl;
1205 filcout << "MUONTRKda : Nb of events used = " << gNEvents << endl;
1210 cout << "\nMUONTRKda : Output logfile generated : " << logOutputFile << endl;
1212 if (gCommand.CompareTo("ped") == 0)
1214 if (!(crocusConfigFile.IsNull()))
1215 cout << "MUONTRKda : CROCUS command file generated : " << crocusOutputFile.Data() << endl;
1217 cout << "MUONTRKda : WARNING no CROCUS command file generated" << endl;
1218 cout << "MUONTRKda : Pedestal Histo file : " << gHistoFileName << endl;
1219 cout << "MUONTRKda : Flat pedestal file (to SHUTTLE) : " << flatOutputFile << endl;
1223 cout << "MUONTRKda : Data file for gain calculation : " << gHistoFileName_gain << endl;
1226 if (gCommand.CompareTo("comp") == 0)
1228 cout << "MUONTRKda : Root Histo. file generated : " << gRootFileName << endl;
1229 cout << "MUONTRKda : Flat gain file (to SHUTTLE) : " << flatOutputFile << endl;
1236 if (gCommand.CompareTo("comp") == 0 || gCommand.CompareTo("ped") == 0)
1238 printf("\n ***** STORE FILE in FES ****** \n");
1240 // to be sure that env variable is set
1241 // gSystem->Setenv("DAQDALIB_PATH", "$DATE_SITE/infoLogger");
1243 if (!flatOutputFile.IsNull())
1246 // flatOutputFile.Prepend("./");
1247 if (gCommand.CompareTo("ped") == 0)
1248 status = daqDA_FES_storeFile(flatOutputFile.Data(),"PEDESTALS");
1250 status = daqDA_FES_storeFile(flatOutputFile.Data(),"GAINS");
1254 printf(" Failed to export file : %d\n",status);
1256 else if(gPrintLevel) printf("Export file: %s\n",flatOutputFile.Data());
1261 printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());