FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairTestDetectorTimeDigiTask.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
9 
10 #include "FairLink.h" // for FairLink
11 #include "FairLogger.h"
12 #include "FairRootManager.h" // for FairRootManager
13 #include "FairTestDetectorDigi.h" // for FairTestDetectorDigi
15 #include "FairTestDetectorPoint.h" // for FairTestDetectorPoint
16 
17 #include <TClonesArray.h> // for TClonesArray
18 #include <TMath.h> // for Sqrt
19 #include <TRandom.h> // for TRandom, gRandom
20 #include <TString.h> // for TString
21 
23  : FairTask()
24  , fTimeResolution(100.)
25  , fPointArray(nullptr)
26  , fDigiArray(nullptr)
27  , fDataBuffer(nullptr)
28  , fTimeOrderedDigi(kFALSE)
29 {}
30 
32 
34 {
36  if (!ioman) {
37  LOG(error) << "FairTestDetectorTimeDigiTask::Init: RootManager not instantiated!";
38  return kFATAL;
39  }
40 
41  fPointArray = static_cast<TClonesArray*>(ioman->GetObject("FairTestDetectorPoint"));
42  if (!fPointArray) {
43  LOG(warn) << "FairTestDetectorTimeDigiTask::Init: No Point array!";
44  return kERROR;
45  }
46 
47  // Create and register output array
48  fDataBuffer = new FairTestDetectorDigiWriteoutBuffer("FairTestDetectorDigi", "TOY", kTRUE);
49  fDataBuffer = static_cast<FairTestDetectorDigiWriteoutBuffer*>(
50  ioman->RegisterWriteoutBuffer("FairTestDetectorDigi", fDataBuffer));
51  fDataBuffer->ActivateBuffering(fTimeOrderedDigi);
52 
53  return kSUCCESS;
54 }
55 
56 void FairTestDetectorTimeDigiTask::Exec(Option_t* /*opt*/)
57 {
58  // fDigiArray->Delete();
59 
60  // fill the map
61  LOG(info) << "EventTime: " << FairRootManager::Instance()->GetEventTime();
62 
63  for (int ipnt = 0; ipnt < fPointArray->GetEntries(); ipnt++) {
64  FairTestDetectorPoint* point = static_cast<FairTestDetectorPoint*>(fPointArray->At(ipnt));
65  if (!point) {
66  continue;
67  }
68 
69  Int_t xPad = CalcPad(point->GetX(), point->GetXOut());
70  Int_t yPad = CalcPad(point->GetY(), point->GetYOut());
71  Int_t zPad = CalcPad(point->GetZ(), point->GetZOut());
72 
73  Double_t timestamp = CalcTimeStamp(point->GetTime());
74 
75  FairTestDetectorDigi* digi = new FairTestDetectorDigi(xPad, yPad, zPad, timestamp);
76  if (fTimeResolution > 0) {
77  digi->SetTimeStampError(fTimeResolution / TMath::Sqrt(fTimeResolution));
78  } else {
79  digi->SetTimeStampError(0);
80  }
81 
82  digi->SetLink(FairLink(-1, FairRootManager::Instance()->GetEntryNr(), "FairTestDetectorPoint", ipnt));
83 
84  Double_t timeOfFlight = point->GetTime();
85  Double_t eventTime = FairRootManager::Instance()->GetEventTime();
86 
87  fDataBuffer->FillNewData(digi, timeOfFlight + eventTime, digi->GetTimeStamp() + 10);
88  }
89 }
90 
91 Int_t FairTestDetectorTimeDigiTask::CalcPad(Double_t posIn, Double_t posOut)
92 {
93  Int_t result = static_cast<Int_t>(posIn + posOut) / 2;
94  return result;
95 }
96 
97 Double_t FairTestDetectorTimeDigiTask::CalcTimeStamp(Double_t timeOfFlight)
98 {
99  Double_t eventTime = FairRootManager::Instance()->GetEventTime();
100  Double_t detectionTime = gRandom->Gaus(0, fTimeResolution);
101 
102  Double_t result = eventTime + timeOfFlight + detectionTime;
103 
104  if (result < 0) {
105  return 0;
106  } else {
107  return result;
108  }
109 }
110 
InitStatus
Definition: FairTask.h:33
Double_t GetZ() const
Definition: FairMCPoint.h:69
Double_t GetX() const
Definition: FairMCPoint.h:67
FairWriteoutBuffer * RegisterWriteoutBuffer(TString branchName, FairWriteoutBuffer *buffer)
virtual void FillNewData(FairTimeStamp *data, double startTime, double activeTime)
static FairRootManager * Instance()
ClassImp(FairEventBuilder)
TObject * GetObject(const char *BrName)
Double_t GetEventTime()
virtual void ActivateBuffering(Bool_t val=kTRUE)
Double_t GetY() const
Definition: FairMCPoint.h:68
Double_t GetTime() const
Definition: FairMCPoint.h:62