/************************************************************************** * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* $Id$ */ // The class which simulates the pulse shape from the PHOS FEE shaper, // make sampled amplitudes, digitize them. // Use case: // AliPHOSPulseGenerator *pulse = new AliPHOSPulseGenerator(energy,time); // Int_t *adcHG = new Int_t[pulse->GetRawFormatTimeBins()]; // Int_t *adcLG= new Int_t[pulse->GetRawFormatTimeBins()]; // pulse->AddNoise(1.); // pulse->MakeSamples(); // pulse->GetSamples(adcHG, adcHG) ; // pulse->Print(); // pulse->Draw(); // // Author: Yuri Kharlov, after Yves Schutz and Per Thomas Hille // --- ROOT system --- #include #include #include #include #include #include // --- AliRoot header files --- #include "AliLog.h" #include "AliPHOSPulseGenerator.h" // --- Standard library --- #include #include using std::cout; using std::endl; ClassImp(AliPHOSPulseGenerator) Double_t AliPHOSPulseGenerator::fgCapa = 1.1; // 1pF Int_t AliPHOSPulseGenerator::fgOrder = 2 ; // order of the Gamma function Double_t AliPHOSPulseGenerator::fgTimePeak = 2.1E-6 ; // tau=2.1 micro seconds Double_t AliPHOSPulseGenerator::fgTimeTrigger = 100E-9 ; // one tick 100 ns Double_t AliPHOSPulseGenerator::fgHighCharge = 8.8; // adjusted for a high gain range of 5 GeV (10 bits) Double_t AliPHOSPulseGenerator::fgHighGain = 6.85; Double_t AliPHOSPulseGenerator::fgHighLowGainFactor = 16.; // adjusted for a low gain range of 80 GeV (10 bits) //----------------------------------------------------------------------------- AliPHOSPulseGenerator::AliPHOSPulseGenerator(Double_t a, Double_t t0) : TObject(), fAmplitude(a), fTZero(t0), fDataHG(0), fDataLG(0), fDigitize(kTRUE) { // Contruct a pulsegenrator object and initializes all necessary parameters // @param a digit amplitude in GeV // @param t0 time delay in nanoseconds of signal relative the first sample. // This value should be between 0 and Ts, where Ts is the sample interval fDataHG = new Double_t[fkTimeBins]; fDataLG = new Double_t[fkTimeBins]; Reset(); } //----------------------------------------------------------------------------- AliPHOSPulseGenerator::AliPHOSPulseGenerator(const AliPHOSPulseGenerator & pulse) : TObject(), fAmplitude(pulse.fAmplitude), fTZero(pulse.fTZero), fDataHG(0), fDataLG(0), fDigitize(kTRUE) { fDataHG = new Double_t[pulse.fkTimeBins]; fDataLG = new Double_t[pulse.fkTimeBins]; for (Int_t i=0; iGaus(0., sigma) ; fDataLG[i] = gRandom->Gaus(0., sigma) ; } } //----------------------------------------------------------------------------- void AliPHOSPulseGenerator::AddNoise(Double_t * /* sigma */, Double_t /* cutoff */) { //Adds correlated Gaussian noise with cutof frequency "cutoff" // @param sigma noise amplitude in entities of ADC levels // @param -30DB cutoff frequency of the noise in entities of sampling frequency AliError("not implemented yet"); } //----------------------------------------------------------------------------- void AliPHOSPulseGenerator::AddPretriggerSamples(Int_t nPresamples) { // Adds pretrigger samples to the sample arrays and replace them // with concatinated and truncated arrays Double_t *tmpDataHG = new Double_t[fkTimeBins]; Double_t *tmpDataLG = new Double_t[fkTimeBins]; Int_t i; for (i=0; i GetRawFormatTimeMax()) signal = 0. ; else { Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder)/fgCapa ; signal = fac * par[2] * TMath::Power(xx/fgTimePeak, fgOrder) * TMath::Exp(-fgOrder*(xx/fgTimePeak)) ; } return signal ; } //__________________________________________________________________ Double_t AliPHOSPulseGenerator::RawResponseFunctionMax(Double_t charge, Double_t gain) { // Maximum value of the shaper response function return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) / ( fgCapa * TMath::Exp(fgOrder) ) ); } //__________________________________________________________________ Bool_t AliPHOSPulseGenerator::MakeSamples() { // for a start time fTZero and an amplitude fAmplitude given by digit, // calculates the raw sampled response AliPHOSPulseGenerator::RawResponseFunction const Int_t kRawSignalOverflow = 0x3FF ; // decimal 1023 Bool_t lowGain = kFALSE ; TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4); for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) { signalF.SetParameter(0, fgHighCharge) ; signalF.SetParameter(1, fgHighGain) ; signalF.SetParameter(2, fAmplitude) ; signalF.SetParameter(3, fTZero) ; Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ; Double_t signal = signalF.Eval(time) ; fDataHG[iTime] += signal; if ( static_cast(fDataHG[iTime]+0.5) > kRawSignalOverflow ){ // larger than 10 bits fDataHG[iTime] = kRawSignalOverflow ; lowGain = kTRUE ; } signalF.SetParameter(0, GetRawFormatLowCharge() ) ; signalF.SetParameter(1, GetRawFormatLowGain() ) ; signal = signalF.Eval(time) ; fDataLG[iTime] += signal; if ( static_cast(fDataLG[iTime]+0.5) > kRawSignalOverflow) // larger than 10 bits fDataLG[iTime] = kRawSignalOverflow ; } // Digitize floating point amplitudes to integers if (fDigitize) Digitize(); return lowGain ; } //__________________________________________________________________ void AliPHOSPulseGenerator::GetSamples(Int_t *adcHG, Int_t *adcLG) const { // Return integer sample arrays adcHG and adcLG for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) { adcHG[iTime] = static_cast(fDataHG[iTime]) ; adcLG[iTime] = static_cast(fDataLG[iTime]) ; } } //__________________________________________________________________ void AliPHOSPulseGenerator::Print(Option_t*) const { // Prints sampled amplitudes to stdout Int_t i; cout << "High gain: "; for (i=0; iSetMarkerStyle(20); graphLG->SetMarkerStyle(20); graphHG->SetMarkerSize(0.4); graphLG->SetMarkerSize(0.4); graphHG->SetTitle("High gain samples"); graphLG->SetTitle("Low gain samples"); TCanvas *c1 = new TCanvas("c1","Raw ALTRO samples",10,10,700,500); c1->SetFillColor(0); c1->Divide(2,1); c1->cd(1); gPad->SetLeftMargin(0.15); graphHG->Draw("AP"); graphHG->GetHistogram()->SetTitleOffset(1.6,"Y"); graphHG->GetHistogram()->SetXTitle("time, #musec"); graphHG->GetHistogram()->SetYTitle("Amplitude, ADC counts"); c1->cd(2); gPad->SetLeftMargin(0.15); graphLG->Draw("AP"); graphLG->GetHistogram()->SetTitleOffset(1.6,"Y"); graphLG->GetHistogram()->SetXTitle("time, #musec"); graphLG->GetHistogram()->SetYTitle("Amplitude, ADC counts"); c1->Update(); }