FairRoot/PandaRoot
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PndTimeStructureAnaTask Class Reference

#include <PndTimeStructureAnaTask.h>

Inheritance diagram for PndTimeStructureAnaTask:

Public Member Functions

 PndTimeStructureAnaTask ()
 
 PndTimeStructureAnaTask (const char *name)
 
virtual ~PndTimeStructureAnaTask ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
virtual void FinishEvent ()
 
virtual void FinishTask ()
 
virtual void SetParContainers ()
 
void AddBranchName (TString name)
 

Protected Member Functions

 ClassDef (PndTimeStructureAnaTask, 1)
 

Protected Attributes

std::vector< DataObjectfData
 
std::vector< DataObjectfDataPrim
 
TH1D * fHistoMixedEvents
 
TH1D * fHistoMixedEventsPrim
 
Int_t fEntryNr
 

Detailed Description

Definition at line 109 of file PndTimeStructureAnaTask.h.

Constructor & Destructor Documentation

PndTimeStructureAnaTask::PndTimeStructureAnaTask ( )
inline

Default constructor

Definition at line 114 of file PndTimeStructureAnaTask.h.

References SetVerbose().

114  :
115  FairTask("TimeStructureAnaTask"), fEntryNr(0)
116  {
117  SetVerbose(0);
118  }
drchit SetVerbose(iVerbose)
PndTimeStructureAnaTask::PndTimeStructureAnaTask ( const char *  name)
inline

Named constructor

Definition at line 121 of file PndTimeStructureAnaTask.h.

References SetVerbose().

121  :
122  FairTask(name), fEntryNr(0)
123  {
124  SetVerbose(0);
125  }
TString name
drchit SetVerbose(iVerbose)
virtual PndTimeStructureAnaTask::~PndTimeStructureAnaTask ( )
inlinevirtual

Destructor

Definition at line 128 of file PndTimeStructureAnaTask.h.

128  {
129  }

Member Function Documentation

void PndTimeStructureAnaTask::AddBranchName ( TString  name)
inline

Definition at line 142 of file PndTimeStructureAnaTask.h.

References fData, and fDataPrim.

142  {
143  fData.push_back(DataObject(name));
144  name += "Prim";
145  fDataPrim.push_back(DataObject(name));
146  }
std::vector< DataObject > fData
std::vector< DataObject > fDataPrim
TString name
PndTimeStructureAnaTask::ClassDef ( PndTimeStructureAnaTask  ,
 
)
protected
void PndTimeStructureAnaTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 72 of file PndTimeStructureAnaTask.cxx.

References fData, fDataPrim, fHistoMixedEvents, fHistoMixedEventsPrim, PndMCTrack::GetMotherID(), and i.

73 {
74  std::cout << " ------------- Event " << FairRootManager::Instance()->GetEntryNr() << " ----------------" << std::endl;
75  bool primaryParticle = false;
76  PndMCTrack* mcTrack = 0;
77  for (size_t i = 0; i < fData.size(); i++){
78  for (int j = 0; j < fData[i].fBranch->GetEntries(); j++){
79  FairTimeStamp* data = (FairTimeStamp*)fData[i].fBranch->At(j);
80  FairMultiLinkedData links = data->GetLinksWithType(FairRootManager::Instance()->GetBranchId("MCTrack"));
81 // std::cout << fData[i].fBranchName.Data() << " " << j << " : " << *data << std::endl;
82  Int_t currentEventNr = -1;
83  if (links.GetNLinks() > 0){
84  currentEventNr = links.GetLink(0).GetEntry();
85  mcTrack = (PndMCTrack*)FairRootManager::Instance()->GetCloneOfLinkData(links.GetLink(0));
86  if (mcTrack != 0)
87  if (mcTrack->GetMotherID() < 0){
88  primaryParticle = true;
89 // std::cout << "Primary Particle found!" << std::endl;
90  }
91  }
92  if (currentEventNr < 0) continue;
93 // if(links.GetNLinks() > 1)
94 // std::cout << fData[i].fBranchName.Data() << " : " << links << std::endl;
95 // std::cout << "CurrentEventNr: " << currentEventNr << " : OldEventNr " << fData[i].fOldEventNr << std::endl;
96  if (fData[i].fEventStartStopMap.count(currentEventNr) == 0){
97  fData[i].fEventStartStopMap[currentEventNr] = std::pair<Double_t, Double_t>(data->GetTimeStamp(), data->GetTimeStamp());
98  }
99  else{
100  fData[i].fEventStartStopMap[currentEventNr].second = data->GetTimeStamp();
101  }
102  if (currentEventNr < fData[i].fOldEventNr){
103  if (fData[i].fEventMixture == kFALSE){
104  fHistoMixedEvents->Fill(i);
105  fData[i].fEventMixture = kTRUE;
106  std::cout << "Event Mixture!" << std::endl;
107  }
108  }
109  else if (currentEventNr == fData[i].fOldEventNr){
110  fData[i].fEventHisto->Fill(data->GetTimeStamp() - fData[i].fOldTimeStamp);
111  fData[i].fOldTimeStamp = data->GetTimeStamp();
112 
113  } else {
114 // std::cout << "HistoFilled: " << data->GetTimeStamp() << " - " << fData[i].fOldTimeStamp << data->GetTimeStamp() - fData[i].fOldTimeStamp << std::endl;
115  fData[i].fTimeHisto->Fill(data->GetTimeStamp() - fData[i].fOldTimeStamp);
116  fData[i].fOldTimeStamp = data->GetTimeStamp();
117  fData[i].fOldEventNr = currentEventNr;
118  fData[i].fEventMixture = kFALSE;
119 
120  }
121 
122  if (primaryParticle){
123  if (fDataPrim[i].fEventStartStopMap.count(currentEventNr) == 0){
124  fDataPrim[i].fEventStartStopMap[currentEventNr] = std::pair<Double_t, Double_t>(data->GetTimeStamp(), data->GetTimeStamp());
125  }
126  else{
127  fDataPrim[i].fEventStartStopMap[currentEventNr].second = data->GetTimeStamp();
128  }
129  if (currentEventNr < fDataPrim[i].fOldEventNr){
130  if (fDataPrim[i].fEventMixture == kFALSE){
131  fHistoMixedEventsPrim->Fill(i);
132  fDataPrim[i].fEventMixture = kTRUE;
133  std::cout << "Event Mixture!" << std::endl;
134  }
135  }
136  else if (currentEventNr == fDataPrim[i].fOldEventNr){
137  fDataPrim[i].fEventHisto->Fill(data->GetTimeStamp() - fDataPrim[i].fOldTimeStamp);
138  fDataPrim[i].fOldTimeStamp = data->GetTimeStamp();
139 
140  } else {
141 // std::cout << "PrimHistoFilled: " << data->GetTimeStamp() << " - " << fDataPrim[i].fOldTimeStamp << data->GetTimeStamp() - fDataPrim[i].fOldTimeStamp << std::endl;
142  fDataPrim[i].fTimeHisto->Fill(data->GetTimeStamp() - fDataPrim[i].fOldTimeStamp);
143  fDataPrim[i].fOldTimeStamp = data->GetTimeStamp();
144  fDataPrim[i].fOldEventNr = currentEventNr;
145  fDataPrim[i].fEventMixture = kFALSE;
146 
147  }
148  }
149  }
150  }
151 }
std::vector< DataObject > fData
Int_t i
Definition: run_full.C:25
std::vector< DataObject > fDataPrim
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
void PndTimeStructureAnaTask::FinishEvent ( )
virtual

Definition at line 155 of file PndTimeStructureAnaTask.cxx.

156 {
157  // fOutputArray->Delete();
158 }
void PndTimeStructureAnaTask::FinishTask ( )
virtual

Definition at line 160 of file PndTimeStructureAnaTask.cxx.

References fData, fDataPrim, fHistoMixedEvents, fHistoMixedEventsPrim, i, and map.

161 {
162  for (size_t i = 0; i < fData.size(); i++){
163  fData[i].fTimeHisto->Write();
164  fData[i].fEventHisto->Write();
165  for (std::map<Int_t, std::pair<Double_t, Double_t> >::iterator iter = fData[i].fEventStartStopMap.begin(); iter != fData[i].fEventStartStopMap.end(); iter++) {
166  fData[i].fEventDiffHisto->Fill(iter->second.second - iter->second.first);
167  }
168  fData[i].fEventDiffHisto->Write();
169  fData[i].FillGraphs();
170  fData[i].fEventGap->Write();
171  fData[i].fOverlap->Write();
172 
173  fData[i].fEventGapPerc->Write();
174  fData[i].fOverlapPerc->Write();
175 
176  fDataPrim[i].fTimeHisto->Write();
177  fDataPrim[i].fEventHisto->Write();
178  for (std::map<Int_t, std::pair<Double_t, Double_t> >::iterator iter = fDataPrim[i].fEventStartStopMap.begin(); iter != fDataPrim[i].fEventStartStopMap.end(); iter++) {
179  fDataPrim[i].fEventDiffHisto->Fill(iter->second.second - iter->second.first);
180  }
181  fDataPrim[i].fEventDiffHisto->Write();
182  fDataPrim[i].FillGraphs();
183  fDataPrim[i].fEventGap->Write();
184  fDataPrim[i].fOverlap->Write();
185 
186  fDataPrim[i].fEventGapPerc->Write();
187  fDataPrim[i].fOverlapPerc->Write();
188  }
189  fHistoMixedEvents->Write();
190  fHistoMixedEventsPrim->Write();
191 }
std::vector< DataObject > fData
Int_t i
Definition: run_full.C:25
std::vector< DataObject > fDataPrim
PndTransMap * map
Definition: sim_emc_apd.C:99
InitStatus PndTimeStructureAnaTask::Init ( )
virtual

Virtual method Init

Definition at line 36 of file PndTimeStructureAnaTask.cxx.

References fData, fDataPrim, fHistoMixedEvents, fHistoMixedEventsPrim, and i.

37 {
38 
39  FairRootManager* ioman = FairRootManager::Instance();
40  if ( ! ioman ) {
41  std::cout << "-E- PndTimeStructureAnaTaskT::Init: "
42  << "RootManager not instantiated!" << std::endl;
43  return kFATAL;
44  }
45 
46  if(fData.size() == 0){
47 
48  fData.push_back(DataObject("MVDSortedPixelDigis"));
49  fData.push_back(DataObject("MVDSortedStripDigis"));
50  fData.push_back(DataObject("STTSortedHits"));
51  fData.push_back(DataObject("GEMDigiSorted"));
52  fData.push_back(DataObject("SciTSortedHit"));
53  fDataPrim.push_back(DataObject("MVDSortedPixelDigisPrim"));
54  fDataPrim.push_back(DataObject("MVDSortedStripDigisPrim"));
55  fDataPrim.push_back(DataObject("STTSortedHitsPrim"));
56  fDataPrim.push_back(DataObject("GEMDigiSortedPrim"));
57  fDataPrim.push_back(DataObject("SciTSortedHitPrim"));
58  }
59 
60  for (size_t i = 0; i < fData.size(); i++){
61  fData[i].fBranch = (TClonesArray*)ioman->GetObject(fData[i].fBranchName);
62  }
63 
64  fHistoMixedEvents = new TH1D("hMixedEvents","hMixedEvents", fData.size(), 0, fData.size());
65  fHistoMixedEventsPrim = new TH1D("hMixedEventsPrim","hMixedEventsPrim", fDataPrim.size(), 0, fDataPrim.size());
66 
67  return kSUCCESS;
68 }
std::vector< DataObject > fData
Int_t i
Definition: run_full.C:25
std::vector< DataObject > fDataPrim
InitStatus PndTimeStructureAnaTask::ReInit ( )
virtual

Definition at line 30 of file PndTimeStructureAnaTask.cxx.

31 {
32  return kSUCCESS;
33 }
virtual void PndTimeStructureAnaTask::SetParContainers ( )
inlinevirtual

Definition at line 140 of file PndTimeStructureAnaTask.h.

140 {};

Member Data Documentation

std::vector<DataObject> PndTimeStructureAnaTask::fData
protected

Definition at line 151 of file PndTimeStructureAnaTask.h.

Referenced by AddBranchName(), Exec(), FinishTask(), and Init().

std::vector<DataObject> PndTimeStructureAnaTask::fDataPrim
protected

Definition at line 152 of file PndTimeStructureAnaTask.h.

Referenced by AddBranchName(), Exec(), FinishTask(), and Init().

Int_t PndTimeStructureAnaTask::fEntryNr
protected

Definition at line 156 of file PndTimeStructureAnaTask.h.

TH1D* PndTimeStructureAnaTask::fHistoMixedEvents
protected

Definition at line 153 of file PndTimeStructureAnaTask.h.

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

TH1D* PndTimeStructureAnaTask::fHistoMixedEventsPrim
protected

Definition at line 154 of file PndTimeStructureAnaTask.h.

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


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