FairRoot/PandaRoot
PndEmcWaveformBuffer.cxx
Go to the documentation of this file.
1 /*
2  * PndEmcWaveformBuffer.cxx
3  */
4 
5 #include "PndEmcWaveformBuffer.h"
6 #include "PndEmcWaveformData.h"
7 #include "PndEmcWaveform.h"
9 #include "PndEmcHit.h"
10 #include "FairLink.h"
11 
12 #include "TClonesArray.h"
13 #include <vector>
14 
15 
16 
18 
19 PndEmcWaveformBuffer::PndEmcWaveformBuffer():FairWriteoutBuffer(), fStoreWaveformData(kFALSE), fWfDataArray(NULL) {
20 }
21 
22 
23 PndEmcWaveformBuffer::PndEmcWaveformBuffer(TString branchName, TString className, TString folderName, Bool_t persistance): FairWriteoutBuffer(branchName, className, folderName, persistance), fStoreWaveformData(kFALSE), fWfDataArray(NULL) {
24 }
25 
26 
28 }
29 
31 
32  FairRootManager* ioman = FairRootManager::Instance();
33 
34  //calculate start and active time for timebased simulation framework
35  //Only times greater or equal than the current EventTime are accepted by the framework, but the absolute time of the first sample might been set to lower times. To allow this, the timebased simulation time parameters are increased as needed.
36  Double_t startTime, activeTime;
37  wfData->GetWaveformSimulator()->GetAbsoluteTimeInterval(wfData, startTime, activeTime);
38  wfData->SetTimeStamp(startTime);
39  wfData->SetTimeOfLastSample(activeTime);
40 
41 // std::cout << "-I- PndEmcWaveformBuffer::FillNewData startTime: " << startTime << " TimeBeforeFirstHit " <<
42 // wfData->GetWaveformSimulator()->GetTimeBeforeFirstHit(wfData) << std::endl;
43 
44  startTime = startTime<ioman->GetEventTime() ? ioman->GetEventTime() : startTime; //shifting startTime towards greater times
45  activeTime += wfData->GetWaveformSimulator()->GetTimeBeforeFirstHit(wfData); // maximal shift in previous step is wfSimulator->GetTimeBeforeFirstHit(), avoid overlapping of generated waves in absolute time domain
46 
47  FairWriteoutBuffer::FillNewData(wfData, startTime, activeTime);
48 
49 }
50 
51 
52 void PndEmcWaveformBuffer::StoreWaveformData(TString branchName, TString folderName, bool persistance) {
53  FairRootManager* ioman = FairRootManager::Instance();
54  fStoreWaveformData = kTRUE;
55  fWfDataBranchName = branchName;
56  ioman->Register(branchName, "PndEmcWaveformData", folderName, persistance);
57  fWfDataArray = ioman->GetTClonesArray(branchName);
58 }
59 
60 
61 
63 
64  FairRootManager* ioman = FairRootManager::Instance();
65  TClonesArray* myArray = ioman->GetTClonesArray(fBranchName);
66 
67  PndEmcWaveformData* wfData = dynamic_cast<PndEmcWaveformData*>(data);
69  wfData->SetEntryNr(FairLink(-1,ioman->GetEntryNr(), fWfDataBranchName, fWfDataArray->GetEntries()));
70  PndEmcWaveform* wave = wfData->GetWaveformSimulator()->Simulate(wfData, myArray);
71 
72  if (fVerbose > 1) {
73  if(wave) {
74  std::cout << "Data Inserted: " << *wave << std::endl;
75  } else {
76  std::cout << "-E in PndEmcWaveformBuffer::AddNewDataToTClonesArray" <<std::endl;
77  }
78  }
79 
80  if(fStoreWaveformData) {
81  new((*fWfDataArray)[fWfDataArray->GetEntries()]) PndEmcWaveformData(*wfData);
82  }
83 }
84 
85 
86 std::vector<std::pair<double, FairTimeStamp*> > PndEmcWaveformBuffer::Modify(std::pair<double, FairTimeStamp*> oldData, std::pair<double, FairTimeStamp* > newData) {
87 
88  PndEmcWaveformData* oldWfData = dynamic_cast<PndEmcWaveformData*>(oldData.second);
89  PndEmcWaveformData* newWfData = dynamic_cast<PndEmcWaveformData*>(newData.second);
90 
91  (*oldWfData)+=(*newWfData);
92  delete newWfData;
93 
94  Double_t startTime, activeTime;
95  oldWfData->GetWaveformSimulator()->GetAbsoluteTimeInterval(oldWfData, startTime, activeTime);
96  oldWfData->SetTimeStamp(startTime);
97  oldWfData->SetTimeOfLastSample(activeTime);
98 
99 
100  activeTime += oldWfData->GetWaveformSimulator()->GetTimeBeforeFirstHit(oldWfData);
101 
102  return std::vector<std::pair<double, FairTimeStamp*> >(1, std::pair<double, FairTimeStamp*>(activeTime, oldWfData));
103 }
104 
105 
106 double PndEmcWaveformBuffer::FindTimeForData(FairTimeStamp* data) {
107  std::map<PndEmcWaveformData, double>::iterator it;
108  PndEmcWaveformData myData = *(PndEmcWaveformData*)data;
109  it = fData_map.find(myData);
110  if (it == fData_map.end())
111  return -1;
112  else
113  return it->second;
114 }
115 
116 void PndEmcWaveformBuffer::FillDataMap(FairTimeStamp* data, double activeTime) {
117  PndEmcWaveformData myData = *(PndEmcWaveformData*)data;
118  fData_map[myData] = activeTime;
119 }
120 
121 
123  FairWriteoutBuffer::DeleteOldData();
124  if(fStoreWaveformData) {
125  fWfDataArray->Delete();
126  }
127 }
128 
129 
131  PndEmcWaveformData myData = *(PndEmcWaveformData*)data;
132  if (fData_map.find(myData) != fData_map.end())
133  fData_map.erase(fData_map.find(myData));
134 }
virtual void EraseDataFromDataMap(FairTimeStamp *data)
std::map< PndEmcWaveformData, double > fData_map
ClassImp(PndEmcWaveformBuffer)
int fVerbose
Definition: poormantracks.C:24
void SetTimeOfLastSample(Double_t time)
virtual Double_t GetTimeBeforeFirstHit(PndEmcWaveformData *)
virtual void StoreWaveformData(TString branchName, TString folderName, bool persistance)
std::vector< std::pair< double, FairTimeStamp * > > Modify(std::pair< double, FairTimeStamp * > oldData, std::pair< double, FairTimeStamp * > newData)
virtual PndEmcWaveform * Simulate(PndEmcWaveformData *wfData, TClonesArray *arrayToStore=NULL)
virtual void AddNewDataToTClonesArray(FairTimeStamp *)
Double_t
PndEmcAbsWaveformSimulator * GetWaveformSimulator()
virtual void GetAbsoluteTimeInterval(PndEmcWaveformData *wfData, Double_t &startTime, Double_t &activeTime)=0
represents a simulated waveform in an emc crystal
virtual double FindTimeForData(FairTimeStamp *data)
buffer for waveforms, used by PndEmcFWEndcapTimebasedWaveforms
virtual void FillNewData(PndEmcWaveformData *)
represents a simulated waveform in an emc crystal, used by PndEmcFWEndcapTimebasedWaveforms ...
virtual void FillDataMap(FairTimeStamp *data, double activeTime)