FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndEmcMultiWaveformToCalibratedDigi Class Reference

Module to take the hit list for the calorimeter and make ADC waveforms from them. More...

#include <PndEmcMultiWaveformToCalibratedDigi.h>

Inheritance diagram for PndEmcMultiWaveformToCalibratedDigi:
PndPersistencyTask

Public Member Functions

 PndEmcMultiWaveformToCalibratedDigi (Int_t verbose=0, Bool_t storedigis=kTRUE)
 
virtual ~PndEmcMultiWaveformToCalibratedDigi ()
 
virtual InitStatus Init ()
 
void SetCalibrationFile (const char *calibrationfilename)
 
virtual void Exec (Option_t *opt)
 
void SetStorageOfData (Bool_t val)
 
virtual void SetPSAAlgorithm (PndEmcAbsPSA *psa)
 
virtual void SetPSAAlgorithmPMT (PndEmcAbsPSA *psa)
 
void SetPersistency (Bool_t val=kTRUE)
 
Bool_t GetPersistency ()
 

Private Member Functions

void ReadCalibrationFile ()
 
virtual void SetParContainers ()
 
 ClassDef (PndEmcMultiWaveformToCalibratedDigi, 2)
 

Private Attributes

TClonesArray * fWaveformArray
 
TClonesArray * fDigiArray
 
Double_t fSampleRate
 
Double_t fSampleRate_PMT
 
Double_t fEnergyDigiThreshold
 
Double_t fASIC_Shaping_int_time
 
Double_t fPMT_Shaping_int_time
 
Double_t fPMT_Shaping_diff_time
 
Double_t fCrystal_time_constant
 
Double_t fShashlyk_time_constant
 
Int_t fNumber_of_samples_in_waveform
 
Int_t fNumber_of_samples_in_waveform_pmt
 
TString fDigiPosMethod
 
Double_t fEmcDigiRescaleFactor
 
Double_t fEmcDigiPositionDepthPWO
 
Double_t fEmcDigiPositionDepthShashlyk
 
PndEmcAbsPulseshapefPulseshape
 
PndEmcAbsPulseshapefPulseshape_pmt
 
PndEmcAbsPSAfpsaAlgorithm
 
PndEmcAbsPSAfpsaAlgorithm_pmt
 
PndEmcDigiParfDigiPar
 
PndEmcRecoParfRecoPar
 
Int_t fVerbose
 
Double_t fWfNormalisation
 
Double_t fWfNormalisation_pmt
 
Double_t fWfNormalisation_proto192
 
std::map< Int_t, Double_tfCalibrationMap
 
std::map< Int_t, Double_tfGainMap
 
TString fCalibrationFileName
 

Detailed Description

Module to take the hit list for the calorimeter and make ADC waveforms from them.

Author
Ch. Hammann chamm.nosp@m.ann@.nosp@m.hiskp.nosp@m..uni.nosp@m.-bonn.nosp@m..de

Definition at line 47 of file PndEmcMultiWaveformToCalibratedDigi.h.

Constructor & Destructor Documentation

PndEmcMultiWaveformToCalibratedDigi::PndEmcMultiWaveformToCalibratedDigi ( Int_t  verbose = 0,
Bool_t  storedigis = kTRUE 
)
PndEmcMultiWaveformToCalibratedDigi::~PndEmcMultiWaveformToCalibratedDigi ( )
virtual

Definition at line 59 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

60 {
61 }

Member Function Documentation

PndEmcMultiWaveformToCalibratedDigi::ClassDef ( PndEmcMultiWaveformToCalibratedDigi  ,
 
)
private
void PndEmcMultiWaveformToCalibratedDigi::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 165 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

References CAMath::Abs(), Bool_t, ctime, Double_t, energy, fCalibrationMap, fDigiArray, fEnergyDigiThreshold, fGainMap, fpsaAlgorithm, fpsaAlgorithm_pmt, fSampleRate, fSampleRate_PMT, fVerbose, fWaveformArray, fWfNormalisation_pmt, fWfNormalisation_proto192, PndEmcWaveform::GetDetectorId(), PndEmcAbsPSA::GetHit(), PndEmcWaveform::GetHitIndex(), PndEmcWaveform::GetModule(), PndEmcMultiWaveform::GetNumberOfWaveforms(), PndEmcWaveform::GetTrackId(), nHits, PndEmcAbsPSA::Process(), rtime, PndEmcMultiWaveform::SetActiveWaveform(), and timer.

166 {
167  TStopwatch timer;
168  if (fVerbose>2){
169  timer.Start();
170  }
171 // Reset output array
172  if ( ! fDigiArray ) Fatal("Exec", "No Digi Array");
173  fDigiArray->Delete();
174  Double_t peakPosition;
176  Double_t digi_time;
177  Int_t i_digi=0; //index of digi in TClonesArray
178  Int_t hitIndex;
179  Int_t detId;
180  Int_t trackId;
181  Int_t module;
182  Int_t nWaveforms = fWaveformArray->GetEntriesFast();
183  Int_t nSignal;
184  Bool_t highgain;
185  PndEmcAbsPSA *thePSA;
186  Double_t theEnergyNorm;
187  Double_t theSampleRate;
188  Int_t nHits[4];
189  Double_t hittimes[4][32];
190  Double_t hitenergies[4][32];
191  std::map<Int_t,Double_t>::iterator it;
192  //cout<<"PndEmcMultiWaveformToCalibratedDigi: "<<nWaveforms<<" waveforms to convert"<<endl;
193  for (Int_t iWaveform=0; iWaveform<nWaveforms; iWaveform++) {
194  PndEmcMultiWaveform* theWaveform = (PndEmcMultiWaveform*) fWaveformArray->At(iWaveform);
195  hitIndex=theWaveform->GetHitIndex();
196  detId=theWaveform->GetDetectorId();
197  trackId=theWaveform->GetTrackId();
198  module=theWaveform->GetModule();
199  nSignal=theWaveform->GetNumberOfWaveforms();
200  if(module == 5){
201  thePSA = fpsaAlgorithm_pmt;
202  theEnergyNorm = fWfNormalisation_pmt;
203  theSampleRate = fSampleRate_PMT;
204  } else {
205  thePSA = fpsaAlgorithm;
206  theEnergyNorm = fWfNormalisation_proto192;;
207  theSampleRate = fSampleRate;
208  }
209  for(Int_t iSignal=0; iSignal < nSignal; iSignal++){
210  theWaveform->SetActiveWaveform(iSignal);
211  nHits[iSignal] = thePSA->Process(theWaveform);
212 // std::cout << "Signal " << iSignal << ": " << nHits[iSignal] << " Hits" << std::endl;
213  for(Int_t iHit = 0; iHit < nHits[iSignal]; iHit ++){
214  thePSA->GetHit(iHit,hitenergies[iSignal][iHit],hittimes[iSignal][iHit]);
215 // std::cout << "Signal " << iSignal << ", Hit " << iHit << ": Energy " << hitenergies[iSignal][iHit] <<std::endl;
216  }
217  }
218 // std::cout << "found "<<nHits[0] << " hits for detid " << detId << std::endl;
219  for(Int_t iHit = 0; iHit < nHits[0]; iHit++){
220  highgain = kTRUE;
221  peakPosition = hittimes[0][iHit];
222  energy = hitenergies[0][iHit];
223  if(nSignal > 1 && energy > 1500.){ //todo make highgain limit a config option
224  //look for hit in lowgain
225 // std::cout << "Hit " << iHit << " is too high (" << energy << ") checking lowgain." << std::endl;
226 // std::cout << "highgain time: " << hittimes[0][iHit] << std::endl;
227  for(Int_t iHit1 = 0; iHit1<nHits[1]; iHit1++){
228 // std::cout << "lowgain hit time: " << hittimes[1][iHit] << std::endl;
229  if(TMath::Abs(peakPosition - hittimes[1][iHit1]) < 10.){ //todo make timediff a config option
230 // std::cout << "Found matching highgain hit " << iHit1 << std::endl;
231  energy = hitenergies[1][iHit1];
232  peakPosition = hittimes[0][iHit];
233  highgain = kFALSE;
234  break;
235  }
236  }
237  }
238 
239 
240  energy/=theEnergyNorm;
241  digi_time = peakPosition/theSampleRate;
242 // std::cout << "looking for calibration of detid " << detId << std::endl;
243  it = fCalibrationMap.find(detId);
244  if(it!=fCalibrationMap.end()){
245 // if(hitIndex == 45) std::cout << "found calibration value of " << it->second << " for hitIndex " << hitIndex << std::endl;
246  energy*=it->second;
247  }else{
248 // if(hitIndex == 45) std::cout << "no calibration value for hitIndex " << hitIndex << std::endl;
249  }
250  if(!highgain){
251 // if(hitIndex == 45) std::cout << "highgain" << std::endl;
252  it=fGainMap.find(detId);
253  if(it!=fGainMap.end()){
254  if(hitIndex == 45) std::cout << "found highgain value of " << it->second << " for hitIndex " << hitIndex << std::endl;
255  energy*=it->second;
256  }
257  }
258  // if(energy>1000){
259  // if(hitIndex == 45) std::cout << "energy: " << energy << " threshold: "<< fEnergyDigiThreshold << endl;
260  // if(hitIndex == 45) std::cout << "creating digi for detid: " << detId << "hitIndex: " << hitIndex <<std::endl;
261  // }
262  if (energy>fEnergyDigiThreshold)
263  {
264 // if(hitIndex == 45) std::cout << "energy: " << energy << endl;
265  PndEmcDigi* myDigi = new((*fDigiArray)[i_digi]) PndEmcDigi(trackId,detId, energy, digi_time, hitIndex);
266  myDigi->AddLink(FairLink("EmcMultiWaveform", iWaveform));
267  i_digi++;
268 
269  }
270  }
271  }
272  if (fVerbose>2){
273  timer.Stop();
274  Double_t rtime = timer.RealTime();
275  Double_t ctime = timer.CpuTime();
276  cout << "PndEmcMultiWaveformToCalibratedDigi, Real time " << rtime << " s, CPU time " << ctime << " s" << endl;
277  }
278 
279 }
Int_t GetNumberOfWaveforms() const
represents the reconstructed hit of one emc crystal
Definition: PndEmcDigi.h:40
virtual Int_t Process(const PndEmcWaveform *waveform)=0
Find Hits in Waveform.
Int_t GetHitIndex() const
virtual void GetHit(Int_t i, Double_t &energy, Double_t &time)=0
Get energy and time of hit.
static T Abs(const T &x)
Definition: PndCAMath.h:39
int nHits
Definition: RiemannTest.C:16
Double_t
int GetTrackId() const
TStopwatch timer
Definition: hit_dirc.C:51
long GetDetectorId() const
Class to hold multiple waveforms from one Emc Hit / ADC readout.
Double_t ctime
Definition: hit_dirc.C:114
Baseclass for pulseshapeanalysis ( featureextraction )
Definition: PndEmcAbsPSA.h:21
Short_t GetModule() const
void SetActiveWaveform(Int_t active=1)
Double_t rtime
Definition: hit_dirc.C:113
Double_t energy
Definition: plot_dirc.C:15
Bool_t PndPersistencyTask::GetPersistency ( )
inlineinherited

Definition at line 32 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency.

Referenced by PndLmdPixelHitProducerFast::GetPersistance(), PndMdtDigitization::Init(), PndMdtHitProducerIdeal::Init(), PndMdtClusterTask::Init(), PndFtsHitProducerRealFast::Init(), PndRichHitProducer::Init(), PndSttHitProducerRealFast::Init(), PndDiscTaskReconstruction::Init(), PndSttHelixHitProducer::Init(), PndDiscTaskPID::Init(), PndIdealTrackFinder::Init(), PndSttMvdGemTracking::Init(), PndMdtTrkProducer::Init(), PndFtsHitProducerRealFull::Init(), PndLmdPixelClusterTask::Init(), PndSttHitProducerRealFull::Init(), PndLmdStripClusterTask::Init(), PndEmcApdHitProducer::Init(), PndMissingPzCleanerTask::Init(), PndEmcMakeRecoHit::Init(), PndEmcMakeClusterOnline::Init(), PndTrackSmearTask::Init(), PndEmcFWEndcapTimebasedWaveforms::Init(), PndSttHitProducerIdeal::Init(), PndEmcFWEndcapDigi::Init(), PndFtsHitProducerIdeal::Init(), PndEmcMakeCluster::Init(), PndMdtPointsToWaveform::Init(), PndDiscTaskDigitization::Init(), PndEmcMakeDigi::Init(), PndSdsTimeWalkCorrTask::Init(), PndLmdPixelHitProducerFast::Init(), PndDrcHitFinder::Init(), PndRichHitFinder::Init(), PndEmcMakeCorr::Init(), PndFtofHitProducerIdeal::Init(), PndEmcHitsToWaveform::Init(), PndSciTDigiTask::Init(), PndDrcHitProducerIdeal::Init(), PndSdsHitProducerIdeal::Init(), PndSciTHitProducerIdeal::Init(), PndEmcHitProducer::Init(), PndRecoMultiKalmanTask2::Init(), PndDrcHitProducerReal::Init(), PndDskFLGHitProducerIdeal::Init(), PndEmcTmpWaveformToDigi::Init(), PndDrcDigiTask::Init(), PndEmcWaveformToDigi::Init(), PndSttMatchTracks::Init(), PndEmcWaveformToCalibratedDigi::Init(), PndTrkTracking2::Init(), PndSttFindTracks::Init(), Init(), PndDrcTimeDigiTask::Init(), PndRecoKalmanTask2::Init(), PndEmcExpClusterSplitter::Init(), PndSdsNoiseProducer::Init(), PndFtsHoughTrackerTask::Init(), PndEmcPhiBumpSplitter::Init(), PndSdsHybridHitProducer::Init(), PndSdsIdealRecoTask::Init(), PndRecoMultiKalmanTask::Init(), PndSdsIdealClusterTask::Init(), PndRecoKalmanTask::Init(), PndSdsStripHitProducerDif::Init(), PndSdsStripHitProducer::Init(), PndGemDigitize::Init(), PndGemFindHits::Init(), PndSdsPixelClusterTask::Init(), PndSdsStripClusterTask::Init(), PndMvdGemTrackFinderOnHits::Init(), PndBarrelTrackFinder::Init(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcMakeBump::PndEmcMakeBump(), PndUnassignedHitsTask::RegisterBranches(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndEmcMakeBump::SetStorageOfData(), and PndEmcFullDigiTask::StoreDigi().

32 { return fPersistency; }
InitStatus PndEmcMultiWaveformToCalibratedDigi::Init ( )
virtual

Virtual method Init

Definition at line 64 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

References PndEmcDigi::depth, Double_t, fASIC_Shaping_int_time, fCalibrationFileName, fCalibrationMap, fCrystal_time_constant, fDigiArray, fDigiPar, fDigiPosMethod, fEmcDigiPositionDepthPWO, fEmcDigiPositionDepthShashlyk, fEmcDigiRescaleFactor, fEnergyDigiThreshold, fNumber_of_samples_in_waveform, fNumber_of_samples_in_waveform_pmt, fPMT_Shaping_diff_time, fPMT_Shaping_int_time, fpsaAlgorithm, fpsaAlgorithm_pmt, fPulseshape, fPulseshape_pmt, fRecoPar, fSampleRate, fSampleRate_PMT, fShashlyk_time_constant, fWaveformArray, fWfNormalisation, fWfNormalisation_pmt, fWfNormalisation_proto192, PndEmcDigiPar::GetASIC_Shaping_int_time(), PndEmcDigiPar::GetCrystal_time_constant(), PndEmcRecoPar::GetEmcDigiPositionDepthPWO(), PndEmcRecoPar::GetEmcDigiPositionDepthShashlyk(), PndEmcDigiPar::GetEnergyDigiThreshold(), PndEmcDigiPar::GetNumber_of_samples_in_waveform(), PndEmcDigiPar::GetNumber_of_samples_in_waveform_pmt(), PndPersistencyTask::GetPersistency(), PndEmcDigiPar::GetPMT_Shaping_diff_time(), PndEmcDigiPar::GetPMT_Shaping_int_time(), PndEmcDigiPar::GetSampleRate(), PndEmcDigiPar::GetSampleRate_PMT(), PndEmcDigiPar::GetShashlyk_time_constant(), PndEmcAbsPSA::Process(), ReadCalibrationFile(), PndEmcDigi::selectDigiPositionMethod(), PndEmcHit::SetEnergy(), PndEmcHit::SetTime(), PndEmcDigi::surface, and PndEmcWaveform::UpdateWaveform().

65 {
66  // Get RootManager
67  FairRootManager* ioman = FairRootManager::Instance();
68  if ( ! ioman )
69  {
70  cout << "-E- PndEmcMultiWaveformToCalibratedDigi::Init: "
71  << "RootManager not instantiated!" << endl;
72  return kFATAL;
73  }
74 
75  // Get input array
76  fWaveformArray = (TClonesArray*) ioman->GetObject("EmcMultiWaveform");
77  if ( ! fWaveformArray ) {
78  cout << "-W- PndEmcMultiWaveformToCalibratedDigi::Init: "
79  << "No PndEmcMultiWaveform array!" << endl;
80  return kERROR;
81  }
82 
83  // Create and register output array
84  fDigiArray = new TClonesArray("PndEmcDigi");
85 
86  ioman->Register("EmcDigi","Emc",fDigiArray,GetPersistency());
99 
100  cout<<"fEmcDigiPositionDepthPWO: "<<fEmcDigiPositionDepthPWO<<endl;
101  cout<<"fEmcDigiPositionDepthShashlyk: "<<fEmcDigiPositionDepthShashlyk<<endl;
102 
103  if (!fDigiPosMethod.CompareTo("surface"))
104  {
106  }
107  else if (!fDigiPosMethod.CompareTo("depth"))
108  {
111  }
112  else
113  {
114  cout << "-W- PndEmcMultiWaveformToCalibratedDigi::Init: "
115  << "Unknown digi position method!" << endl;
116  return kERROR;
117  }
118 
121 
122 
123  if(fpsaAlgorithm == NULL){
124  // Matched digital filter
125  // Parameters of the filter are hardcoded at the moment
126  // For different pulseshape in barrel and endcaps different filters should be implemented
127  std::vector<Double_t> params;
128  params.push_back(30); // width
129  params.push_back(fSampleRate); // Sample rate
131  }
132 
133  if(fpsaAlgorithm_pmt == NULL){
134  // Pulse shape analysis algorithm.
135  // Simple parabolic fit.
137  }
138 
139  // Determine normalisation constant for PndEmcMultiWaveform
140  PndEmcWaveform *tmpwaveform=new PndEmcWaveform(0,101010001, fNumber_of_samples_in_waveform);
141  PndEmcWaveform *tmpwaveform2=new PndEmcWaveform(0,101010001, fNumber_of_samples_in_waveform_pmt);
142 
143  PndEmcHit *gevHit=new PndEmcHit();
144  gevHit->SetEnergy(1.0);
145  gevHit->SetTime(0.);
146  tmpwaveform->UpdateWaveform(gevHit, 0, false, 1., 0., fSampleRate, fPulseshape);
147  tmpwaveform2->UpdateWaveform(gevHit, 0, false, 1., 0., fSampleRate_PMT, fPulseshape_pmt);
148  Double_t tmpPeakPosition;
149  Double_t tmpPeakPosition2;
150  fpsaAlgorithm->Process(tmpwaveform,fWfNormalisation,tmpPeakPosition);
151  fpsaAlgorithm_pmt->Process(tmpwaveform2,fWfNormalisation_pmt,tmpPeakPosition2);
152 
154 
155  if(fCalibrationFileName !="")
157 
158  cout << "-I- PndEmcMultiWaveformToCalibratedDigi: Read "<< fCalibrationMap.size() << " Calibration Entries" << endl;
159 
160  cout << "-I- PndEmcMultiWaveformToCalibratedDigi: Intialization successfull" << endl;
161 
162  return kSUCCESS;
163 }
Double_t GetEmcDigiPositionDepthShashlyk()
Definition: PndEmcRecoPar.h:25
Double_t GetSampleRate()
Definition: PndEmcDigiPar.h:39
Double_t GetPMT_Shaping_diff_time()
Definition: PndEmcDigiPar.h:33
static void selectDigiPositionMethod(PositionMethod, double positionDepthPWO=0., double positionDepthShahslyk=0., double rescaleFactor=1.)
Definition: PndEmcDigi.cxx:153
Int_t GetNumber_of_samples_in_waveform()
Definition: PndEmcDigiPar.h:44
Pulseshape analysis for ADC waveforms.
virtual Int_t Process(const PndEmcWaveform *waveform)=0
Find Hits in Waveform.
Double_t GetCrystal_time_constant()
Definition: PndEmcDigiPar.h:36
Pulseshape analysis for ADC waveforms.
Double_t GetASIC_Shaping_int_time()
Definition: PndEmcDigiPar.h:31
Pulseshape from an CRRC-Shaper.
Double_t GetEnergyDigiThreshold()
Definition: PndEmcDigiPar.h:42
Pulseshape from an APFEL ASIC preamplifier shaper.
Double_t
void UpdateWaveform(PndEmcHit *hit, Double_t pePerMeV, Bool_t usePhotonStatistic, Double_t excessNoiseFactor, Double_t firstADCBinTime, Double_t sampleRate, PndEmcAbsPulseshape *pulseshape, Double_t=0)
virtual void SetEnergy(Double32_t energy)
Definition: PndEmcHit.h:50
represents a simulated waveform in an emc crystal
represents the deposited energy of one emc crystal from simulation
Definition: PndEmcHit.h:26
Double_t GetSampleRate_PMT()
Definition: PndEmcDigiPar.h:40
virtual void SetTime(Double32_t time)
Definition: PndEmcHit.h:51
Double_t GetShashlyk_time_constant()
Definition: PndEmcDigiPar.h:37
Int_t GetNumber_of_samples_in_waveform_pmt()
Definition: PndEmcDigiPar.h:45
Double_t GetPMT_Shaping_int_time()
Definition: PndEmcDigiPar.h:32
Double_t GetEmcDigiPositionDepthPWO()
Definition: PndEmcRecoPar.h:24
void PndEmcMultiWaveformToCalibratedDigi::ReadCalibrationFile ( )
private

Definition at line 304 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

References buf, fCalibrationFileName, fCalibrationMap, fGainMap, fVerbose, and TString.

Referenced by Init().

304  {
305  std::ifstream in(fCalibrationFileName,std::ifstream::in);
306 // in.open(fCalibrationFileName);
307  if (!in.good()) {
308  std::cerr << "Cannot open calibration file!" << endl;
309  return;
310  }
311 
312  char buf[255];
313  while (in.getline(buf, 255)) {
314  TString tmp = buf;
315  TObjArray *tokens = tmp.Tokenize(" \t;");
316  if(tokens->GetEntries()<2){
317  continue;
318  }
319  tmp=tokens->UncheckedAt(0)->GetName();
320  TString tmp2 = tokens->UncheckedAt(1)->GetName();
321  if(tmp.IsDigit() && tmp2.IsFloat()){
322  fCalibrationMap.insert(std::pair<Int_t,Double_t>(tmp.Atoi(),tmp2.Atof()));
323 
324  }
325  if(tokens->GetEntries()<3){
326  continue;
327  }
328  tmp2 = tokens->UncheckedAt(2)->GetName();
329  if(tmp.IsDigit() && tmp2.IsFloat()){
330  fGainMap.insert(std::pair<Int_t,Double_t>(tmp.Atoi(),tmp2.Atof()));
331 
332  }
333  delete tokens;
334 
335 
336  }
337  if(fVerbose >1){
338  std::cout << "calibration values:" << std::endl;
339  std::map<Int_t,Double_t>::iterator it;
340  std::map<Int_t,Double_t>::iterator it2;
341  Int_t detid;
342  for ( it=fCalibrationMap.begin() ; it != fCalibrationMap.end(); it++ ){
343  detid = (*it).first;
344  std::cout << detid << "\t" << (*it).second << "\t";
345  it2=fGainMap.find(detid);
346  if(it2!=fGainMap.end()){
347  std::cout << (*it2).second;
348  }
349  std::cout << endl;
350  }
351  }
352 }
return buf
void PndEmcMultiWaveformToCalibratedDigi::SetCalibrationFile ( const char *  calibrationfilename)
inline
void PndEmcMultiWaveformToCalibratedDigi::SetParContainers ( )
privatevirtual

Reconstruction parameter container Get parameter containers

Definition at line 281 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

References fDigiPar, fRecoPar, and run.

281  {
282 
283  // Get run and runtime database
284  FairRun* run = FairRun::Instance();
285  if ( ! run ) Fatal("SetParContainers", "No analysis run");
286 
287  FairRuntimeDb* db = run->GetRuntimeDb();
288  if ( ! db ) Fatal("SetParContainers", "No runtime database");
289 
290  // Get Emc digitisation parameter container
291  fDigiPar = (PndEmcDigiPar*) db->getContainer("PndEmcDigiPar");
292 
293  // Get Emc reconstruction parameter container
294  fRecoPar = (PndEmcRecoPar*) db->getContainer("PndEmcRecoPar");
295 
296 }
Int_t run
Definition: autocutx.C:47
parameter set of Emc digitisation
Definition: PndEmcDigiPar.h:12
Parameter set for Emc Reco.
Definition: PndEmcRecoPar.h:12
void PndPersistencyTask::SetPersistency ( Bool_t  val = kTRUE)
inlineinherited

Definition at line 31 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency, and val.

Referenced by barrelTrackFinder(), digi_complete(), digi_complete_newSTT(), digiOnly_complete(), PndBarrelTrackFinder::PndBarrelTrackFinder(), PndCATracking::PndCATracking(), PndDrcHitFinder::PndDrcHitFinder(), PndEmc2DLocMaxFinder::PndEmc2DLocMaxFinder(), PndEmcExpClusterSplitter::PndEmcExpClusterSplitter(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcFWEndcapDigi::PndEmcFWEndcapDigi(), PndEmcFWEndcapTimebasedWaveforms::PndEmcFWEndcapTimebasedWaveforms(), PndEmcHitProducer::PndEmcHitProducer(), PndEmcHitsToWaveform::PndEmcHitsToWaveform(), PndEmcMakeBump::PndEmcMakeBump(), PndEmcMakeCluster::PndEmcMakeCluster(), PndEmcMakeClusterOnline::PndEmcMakeClusterOnline(), PndEmcMakeDigi::PndEmcMakeDigi(), PndEmcMakeRecoHit::PndEmcMakeRecoHit(), PndEmcMultiWaveformToCalibratedDigi(), PndEmcPhiBumpSplitter::PndEmcPhiBumpSplitter(), PndEmcTmpWaveformToDigi::PndEmcTmpWaveformToDigi(), PndEmcWaveformToCalibratedDigi::PndEmcWaveformToCalibratedDigi(), PndEmcWaveformToDigi::PndEmcWaveformToDigi(), PndFtofHitProducerIdeal::PndFtofHitProducerIdeal(), PndFtsCATracking::PndFtsCATracking(), PndFtsHitProducerIdeal::PndFtsHitProducerIdeal(), PndFtsHitProducerRealFast::PndFtsHitProducerRealFast(), PndFtsHitProducerRealFull::PndFtsHitProducerRealFull(), PndFtsHoughTrackerTask::PndFtsHoughTrackerTask(), PndGemDigitize::PndGemDigitize(), PndGemFindHits::PndGemFindHits(), PndIdealTrackFinder::PndIdealTrackFinder(), PndLmdPixelClusterTask::PndLmdPixelClusterTask(), PndLmdPixelHitProducerFast::PndLmdPixelHitProducerFast(), PndMdtClusterTask::PndMdtClusterTask(), PndMdtDigitization::PndMdtDigitization(), PndMdtHitProducerIdeal::PndMdtHitProducerIdeal(), PndMdtPointsToWaveform::PndMdtPointsToWaveform(), PndMdtTrkProducer::PndMdtTrkProducer(), PndMissingPzCleanerTask::PndMissingPzCleanerTask(), PndMvdGemTrackFinderOnHits::PndMvdGemTrackFinderOnHits(), PndMvdHitProducerIdeal::PndMvdHitProducerIdeal(), PndMvdPixelClusterTask::PndMvdPixelClusterTask(), PndMvdTimeWalkCorrTask::PndMvdTimeWalkCorrTask(), PndMvdToPix4ClusterTask::PndMvdToPix4ClusterTask(), PndRecoKalmanTask::PndRecoKalmanTask(), PndRecoKalmanTask2::PndRecoKalmanTask2(), PndRecoMultiKalmanTask::PndRecoMultiKalmanTask(), PndRecoMultiKalmanTask2::PndRecoMultiKalmanTask2(), PndRichHitFinder::PndRichHitFinder(), PndRichHitProducer::PndRichHitProducer(), PndSciTDigiTask::PndSciTDigiTask(), PndSciTHitProducerIdeal::PndSciTHitProducerIdeal(), PndSdsHitProducerIdeal::PndSdsHitProducerIdeal(), PndSdsHybridHitProducer::PndSdsHybridHitProducer(), PndSdsIdealClusterTask::PndSdsIdealClusterTask(), PndSdsIdealRecoTask::PndSdsIdealRecoTask(), PndSdsNoiseProducer::PndSdsNoiseProducer(), PndSdsPixelClusterTask::PndSdsPixelClusterTask(), PndSdsStripClusterTask::PndSdsStripClusterTask(), PndSdsStripHitProducer::PndSdsStripHitProducer(), PndSdsTimeWalkCorrTask::PndSdsTimeWalkCorrTask(), PndSttFindTracks::PndSttFindTracks(), PndSttHelixHitProducer::PndSttHelixHitProducer(), PndSttHitProducerIdeal::PndSttHitProducerIdeal(), PndSttHitProducerRealFast::PndSttHitProducerRealFast(), PndSttHitProducerRealFull::PndSttHitProducerRealFull(), PndSttMatchTracks::PndSttMatchTracks(), PndSttMvdGemTracking::PndSttMvdGemTracking(), PndTrackSmearTask::PndTrackSmearTask(), PndTrkTracking2::PndTrkTracking2(), reco(), reco_complete(), reco_complete_gf2(), reco_complete_newSTT(), reco_complete_sec(), recoideal_complete(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndLmdPixelHitProducerFast::SetPersistance(), PndSdsHitProducerIdeal::SetPersistance(), PndSttMvdGemTracking::SetPersistenc(), PndMdtClusterTask::SetPersistence(), PndSttHelixHitProducer::SetPersistence(), PndMissingPzCleanerTask::SetPersistence(), PndFtsHitProducerRealFast::SetPersistence(), PndFtsHitProducerRealFull::SetPersistence(), PndSttHitProducerIdeal::SetPersistence(), PndSttHitProducerRealFull::SetPersistence(), PndSttHitProducerRealFast::SetPersistence(), PndFtsHitProducerIdeal::SetPersistence(), PndTrackSmearTask::SetPersistence(), PndSciTHitProducerIdeal::SetPersistence(), PndIdealTrackFinder::SetPersistence(), PndSttMatchTracks::SetPersistence(), PndSttFindTracks::SetPersistence(), PndFtsHoughTrackerTask::SetPersistence(), PndTrkTracking2::SetPersistence(), PndEmcMakeRecoHit::SetStorageOfData(), PndEmcMakeClusterOnline::SetStorageOfData(), PndEmcFWEndcapDigi::SetStorageOfData(), PndEmcFWEndcapTimebasedWaveforms::SetStorageOfData(), PndEmcMakeDigi::SetStorageOfData(), PndMdtPointsToWaveform::SetStorageOfData(), PndEmc2DLocMaxFinder::SetStorageOfData(), PndEmcMakeCluster::SetStorageOfData(), PndEmcHitsToWaveform::SetStorageOfData(), PndEmcMakeBump::SetStorageOfData(), PndEmcTmpWaveformToDigi::SetStorageOfData(), PndEmcWaveformToDigi::SetStorageOfData(), PndEmcWaveformToCalibratedDigi::SetStorageOfData(), SetStorageOfData(), PndEmcExpClusterSplitter::SetStorageOfData(), PndEmcPhiBumpSplitter::SetStorageOfData(), standard_tracking(), and PndEmcFullDigiTask::StoreDigi().

31 { fPersistency = val; }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
virtual void PndEmcMultiWaveformToCalibratedDigi::SetPSAAlgorithm ( PndEmcAbsPSA psa)
inlinevirtual
virtual void PndEmcMultiWaveformToCalibratedDigi::SetPSAAlgorithmPMT ( PndEmcAbsPSA psa)
inlinevirtual
void PndEmcMultiWaveformToCalibratedDigi::SetStorageOfData ( Bool_t  val)

Definition at line 298 of file PndEmcMultiWaveformToCalibratedDigi.cxx.

References PndPersistencyTask::SetPersistency().

299 {
301  return;
302 }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
void SetPersistency(Bool_t val=kTRUE)

Member Data Documentation

Double_t PndEmcMultiWaveformToCalibratedDigi::fASIC_Shaping_int_time
private

Definition at line 84 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

TString PndEmcMultiWaveformToCalibratedDigi::fCalibrationFileName
private
std::map<Int_t,Double_t> PndEmcMultiWaveformToCalibratedDigi::fCalibrationMap
private

Definition at line 115 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), Init(), and ReadCalibrationFile().

Double_t PndEmcMultiWaveformToCalibratedDigi::fCrystal_time_constant
private

Definition at line 87 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

TClonesArray* PndEmcMultiWaveformToCalibratedDigi::fDigiArray
private

output array of EmcDigis

Definition at line 79 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

PndEmcDigiPar* PndEmcMultiWaveformToCalibratedDigi::fDigiPar
private

Definition at line 103 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init(), and SetParContainers().

TString PndEmcMultiWaveformToCalibratedDigi::fDigiPosMethod
private
Double_t PndEmcMultiWaveformToCalibratedDigi::fEmcDigiPositionDepthPWO
private

Definition at line 95 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fEmcDigiPositionDepthShashlyk
private

Definition at line 96 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fEmcDigiRescaleFactor
private
Double_t PndEmcMultiWaveformToCalibratedDigi::fEnergyDigiThreshold
private

Definition at line 83 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

std::map<Int_t,Double_t> PndEmcMultiWaveformToCalibratedDigi::fGainMap
private

Definition at line 116 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and ReadCalibrationFile().

Int_t PndEmcMultiWaveformToCalibratedDigi::fNumber_of_samples_in_waveform
private

Definition at line 89 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Int_t PndEmcMultiWaveformToCalibratedDigi::fNumber_of_samples_in_waveform_pmt
private

Definition at line 90 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fPMT_Shaping_diff_time
private

Definition at line 86 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fPMT_Shaping_int_time
private

Definition at line 85 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

PndEmcAbsPSA* PndEmcMultiWaveformToCalibratedDigi::fpsaAlgorithm
private

Definition at line 100 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), Init(), and SetPSAAlgorithm().

PndEmcAbsPSA* PndEmcMultiWaveformToCalibratedDigi::fpsaAlgorithm_pmt
private

Definition at line 101 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), Init(), and SetPSAAlgorithmPMT().

PndEmcAbsPulseshape* PndEmcMultiWaveformToCalibratedDigi::fPulseshape
private

Definition at line 98 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

PndEmcAbsPulseshape* PndEmcMultiWaveformToCalibratedDigi::fPulseshape_pmt
private

Definition at line 99 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

PndEmcRecoPar* PndEmcMultiWaveformToCalibratedDigi::fRecoPar
private

Digitisation parameter container

Definition at line 104 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init(), and SetParContainers().

Double_t PndEmcMultiWaveformToCalibratedDigi::fSampleRate
private

Definition at line 81 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fSampleRate_PMT
private

Definition at line 82 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fShashlyk_time_constant
private

Definition at line 88 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Int_t PndEmcMultiWaveformToCalibratedDigi::fVerbose
private

Verbosity level

Definition at line 109 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), PndEmcMultiWaveformToCalibratedDigi(), and ReadCalibrationFile().

TClonesArray* PndEmcMultiWaveformToCalibratedDigi::fWaveformArray
private

Input array of PndEmcMultiWaveforms

Definition at line 76 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fWfNormalisation
private

Definition at line 111 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fWfNormalisation_pmt
private

Definition at line 112 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().

Double_t PndEmcMultiWaveformToCalibratedDigi::fWfNormalisation_proto192
private

Definition at line 113 of file PndEmcMultiWaveformToCalibratedDigi.h.

Referenced by Exec(), and Init().


The documentation for this class was generated from the following files: