FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairTestDetectorRecoTask.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
14 #include "FairTestDetectorHit.h" // for FairTestDetectorHit
15 
16 #include <TClonesArray.h> // for TClonesArray
17 #include <TMath.h> // for Sqrt
18 #include <TVector3.h> // for TVector3
19 
21  : FairTask()
22  , fDigiArray(nullptr)
23  , fHitArray(nullptr)
24 {}
25 
27  : FairTask()
28  , fDigiArray(nullptr)
29  , fHitArray(nullptr)
30 {
31  fVerbose = verbose;
32 }
33 
35 
37 {
39  if (!ioman) {
40  LOG(error) << "-E- FairTestDetectorRecoTask::Init: RootManager not instantiated!";
41  return kFATAL;
42  }
43 
44  fDigiArray = static_cast<TClonesArray*>(ioman->GetObject("FairTestDetectorDigi"));
45  if (!fDigiArray) {
46  LOG(warn) << "-W- FairTestDetectorRecoTask::Init: No Point array!";
47  return kERROR;
48  }
49 
50  // Create and register output array
51  fHitArray = new TClonesArray("FairTestDetectorHit");
52  ioman->Register("FairTestDetectorHit", "FairTestDetector", fHitArray, kTRUE);
53 
54  return kSUCCESS;
55 }
56 
57 void FairTestDetectorRecoTask::Exec(Option_t* /*opt*/)
58 {
59  fHitArray->Clear();
60 
61  // fill the map
62  for (int ipnt = 0; ipnt < fDigiArray->GetEntriesFast(); ipnt++) {
63  FairTestDetectorDigi* digi = static_cast<FairTestDetectorDigi*>(fDigiArray->At(ipnt));
64  if (!digi)
65  continue;
66 
67  /*
68  LOG(debug) << " x= " << digi->GetX()
69  << " y= " << digi->GetY()
70  << " z= " << digi->GetZ()
71  << " t= " << digi->GetTimeStamp();
72  // */
73 
74  TVector3 pos(digi->GetX() + 0.5, digi->GetY() + 0.5, digi->GetZ() + 0.5);
75  TVector3 dpos(1 / TMath::Sqrt(12), 1 / TMath::Sqrt(12), 1 / TMath::Sqrt(12));
76 
77  FairTestDetectorHit* hit = new ((*fHitArray)[ipnt]) FairTestDetectorHit(-1, -1, pos, dpos);
78  if (fStreamProcessing == kFALSE)
79  hit->AddLink(FairLink(-1,
80  FairRootManager::Instance()->GetEntryNr(),
81  FairRootManager::Instance()->GetBranchId("FairTestDetectorDigi"),
82  ipnt));
83  hit->SetTimeStamp(digi->GetTimeStamp());
85  }
86 }
87 
Bool_t fStreamProcessing
Definition: FairTask.h:103
InitStatus
Definition: FairTask.h:33
static FairRootManager * Instance()
ClassImp(FairEventBuilder)
TObject * GetObject(const char *BrName)
Int_t fVerbose
Definition: FairTask.h:100
Double_t GetTimeStampError() const
Definition: FairTimeStamp.h:45
Double_t GetTimeStamp() const
Definition: FairTimeStamp.h:44
void SetTimeStamp(Double_t t)
Definition: FairTimeStamp.h:47
virtual void Exec(Option_t *opt)
void Register(const char *name, const char *Foldername, TNamed *obj, Bool_t toFile)
void SetTimeStampError(Double_t t)
Definition: FairTimeStamp.h:48