CbmRoot
NicaUnigenEventInterface.cxx
Go to the documentation of this file.
1 /*
2  * NicaUnigenSource.cxx
3  *
4  * Created on: 2 sie 2017
5  * Author: Daniel Wielanek
6  * E-mail: daniel.wielanek@gmail.com
7  * Warsaw University of Technology, Faculty of Physics
8  */
10 #include "FairRootManager.h"
12 
14  fEvent = new UEvent();
15 }
16 
17 void NicaUnigenEventInterface::Compress(Int_t* map, Int_t map_size) {
18  Int_t track_pos = 0;
19  for (int i = 0; i < map_size; i++) {
20  Int_t good_track = map[i];
21  for (int j = track_pos; j < good_track; j++) {
22  fEvent->RemoveAt(j);
23  }
24  track_pos = good_track + 1;
25  }
26 }
27 
28 void NicaUnigenEventInterface::CopyData(NicaEventInterface* s) {
29 #ifdef UNIGEN_OLD
30  CopyUnigen(((NicaUnigenEventInterface*) s)->fEvent, fEvent);
31 #else
32  *fEvent = *((NicaUnigenEventInterface*) s)->fEvent;
33 #endif
34 }
35 
36 void NicaUnigenEventInterface::CopyAndCompress(NicaEventInterface* s,
37  Int_t* map,
38  Int_t map_size) {
40  fEvent->SetB(ev->fEvent->GetB());
41  fEvent->SetPhi(ev->fEvent->GetPhi());
42  fEvent->SetNes(ev->fEvent->GetNes());
44  fEvent->SetStepT(ev->fEvent->GetStepT());
45 #ifdef UNIGEN_OLD
46  fEvent->GetParticleList()->Clear();
47 #else
48  TString comment;
49  ev->fEvent->GetComment(comment);
50  fEvent->SetComment(comment);
51  fEvent->Clear();
52 #endif
53  for (int i = 0; i < map_size; i++) {
54  fEvent->AddParticle(*ev->fEvent->GetParticle(map[i]));
55  }
56 }
57 
59  FairRootManager* manager = FairRootManager::Instance();
60  if (CanDeleteEvent()) {
61  if (fEvent) delete fEvent;
62  }
63  fEvent = (UEvent*) manager->GetObject("UEvent.");
64 }
65 
66 void NicaUnigenEventInterface::Boost(Double_t vx, Double_t vy, Double_t vz) {
67  for (int i = 0; i < fEvent->GetNpa(); i++) {
69  TLorentzVector mom = p->GetMomentum();
70  TLorentzVector pos = p->GetPosition();
71  mom.Boost(vx, vy, vz);
72  pos.Boost(vx, vy, vz);
73  p->SetMomentum(mom);
74  p->SetPosition(pos);
75  }
76 }
77 
79  if (CanDeleteEvent()) {
80  if (fEvent) delete fEvent;
81  }
82 }
83 
84 NicaTrackInterface* NicaUnigenEventInterface::GetTrackInterface() const {
85  return new NicaUnigenTrackInterface();
86 }
87 
89  if (fEvent == NULL) fEvent = new UEvent();
90  FairRootManager* manager = FairRootManager::Instance();
91  manager->Register("Event", "", (TNamed*) fEvent, write);
92 }
93 
94 void NicaUnigenEventInterface::FillTrackInterface(NicaTrackInterface* track,
95  Int_t index) {
96  track->SetRawTrack(fEvent->GetParticle(index));
97 }
98 #ifdef UNIGEN_OLD
99 void NicaUnigenEventInterface::CopyUnigen(UEvent* from, UEvent* to) {
100  to->GetParticleList()->Clear();
101  to->SetB(from->GetB());
102  to->SetPhi(from->GetPhi());
103  to->SetNes(from->GetNes());
104  to->SetStepNr(from->GetStepNr());
105  to->SetStepT(from->GetStepT());
106  for (int i = 0; i < from->GetNpa(); i++) {
107  to->AddParticle(*from->GetParticle(i));
108  }
109 }
110 #endif
UEvent::GetParticleList
TClonesArray * GetParticleList() const
Definition: UEvent.h:38
UEvent::GetPhi
Double_t GetPhi() const
Definition: UEvent.h:33
UEvent::GetB
Double_t GetB() const
Definition: UEvent.h:32
UParticle
Definition: UParticle.h:10
NicaUnigenEventInterface::GetTrackInterface
virtual NicaTrackInterface * GetTrackInterface() const
Definition: NicaUnigenEventInterface.cxx:84
NicaUnigenEventInterface::ConnectToTree
virtual void ConnectToTree()
Definition: NicaUnigenEventInterface.cxx:58
UEvent::GetNpa
Int_t GetNpa() const
Definition: UEvent.h:37
NicaUnigenEventInterface::FillTrackInterface
virtual void FillTrackInterface(NicaTrackInterface *track, Int_t index)
Definition: NicaUnigenEventInterface.cxx:94
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
UEvent::GetComment
void GetComment(TString &comment) const
Definition: UEvent.h:39
NicaUnigenEventInterface::CopyData
virtual void CopyData(NicaEventInterface *s)
Definition: NicaUnigenEventInterface.cxx:28
UParticle::SetMomentum
void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: UParticle.h:115
NicaUnigenTrackInterface
Definition: NicaUnigenTrackInterface.h:17
UEvent::AddParticle
void AddParticle(Int_t index, Int_t pdg, Int_t status, Int_t parent, Int_t parentDecay, Int_t mate, Int_t decay, Int_t child[2], Double_t px, Double_t py, Double_t pz, Double_t e, Double_t x, Double_t y, Double_t z, Double_t t, Double_t weight)
Definition: UEvent.cxx:103
UEvent::GetStepNr
Int_t GetStepNr() const
Definition: UEvent.h:35
UEvent
Definition: UEvent.h:12
UEvent::RemoveAt
void RemoveAt(Int_t i)
Definition: UEvent.cxx:230
NicaUnigenEventInterface::CopyAndCompress
virtual void CopyAndCompress(NicaEventInterface *s, Int_t *map, Int_t map_size)
Definition: NicaUnigenEventInterface.cxx:36
UEvent::SetStepNr
void SetStepNr(Int_t stepNr)
Definition: UEvent.h:52
NicaUnigenEventInterface.h
NicaUnigenEventInterface::Boost
virtual void Boost(Double_t vx, Double_t vy, Double_t vz)
Definition: NicaUnigenEventInterface.cxx:66
UEvent::SetB
void SetB(Double_t b)
Definition: UEvent.h:49
UEvent::SetNes
void SetNes(Int_t nes)
Definition: UEvent.h:51
UEvent::SetStepT
void SetStepT(Double_t stepT)
Definition: UEvent.h:53
NicaUnigenEventInterface::~NicaUnigenEventInterface
virtual ~NicaUnigenEventInterface()
Definition: NicaUnigenEventInterface.cxx:78
UEvent::SetComment
void SetComment(const char *comment)
Definition: UEvent.h:54
UEvent::GetParticle
UParticle * GetParticle(Int_t index) const
Definition: UEvent.cxx:92
NicaUnigenEventInterface::NicaUnigenEventInterface
NicaUnigenEventInterface()
Definition: NicaUnigenEventInterface.cxx:13
NicaUnigenTrackInterface.h
NicaUnigenEventInterface::Compress
virtual void Compress(Int_t *map, Int_t map_size)
Definition: NicaUnigenEventInterface.cxx:17
UParticle::GetMomentum
TLorentzVector GetMomentum() const
Definition: UParticle.h:81
UParticle::GetPosition
TLorentzVector GetPosition() const
Definition: UParticle.h:91
UParticle::SetPosition
void SetPosition(Double_t x, Double_t y, Double_t z, Double_t t)
Definition: UParticle.h:131
UEvent::GetStepT
Double_t GetStepT() const
Definition: UEvent.h:36
NicaUnigenEventInterface
Definition: NicaUnigenEventInterface.h:21
pos
TVector3 pos
Definition: CbmMvdSensorDigiToHitTask.cxx:60
NicaUnigenEventInterface::fEvent
UEvent * fEvent
Definition: NicaUnigenEventInterface.h:23
UEvent::GetNes
Int_t GetNes() const
Definition: UEvent.h:34
NicaUnigenEventInterface::Register
virtual void Register(Bool_t write)
Definition: NicaUnigenEventInterface.cxx:88
UEvent::SetPhi
void SetPhi(Double_t phi)
Definition: UEvent.h:50
UEvent::Clear
void Clear(Option_t *="")
Definition: UEvent.cxx:221