CbmRoot
Loading...
Searching...
No Matches
HalCbmUnigenEvent.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023-2023 Warsaw University of Technology, Warsaw
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Daniel Wielanek [committer] */
4#include "HalCbmUnigenEvent.h"
5
7
8#include <TClonesArray.h>
9#include <TDatabasePDG.h>
10#include <TLorentzVector.h>
11#include <TParticlePDG.h>
12
13#include <Hal/DataManager.h>
14#include <Hal/EventInterface.h>
15#include <Hal/McTrack.h>
16#include <Hal/Track.h>
17
18HalCbmUnigenEvent::HalCbmUnigenEvent() : Hal::McEvent("HalCbmUnigenTrack") {}
19
20void HalCbmUnigenEvent::Update(Hal::EventInterface* interface)
21{
22 UEvent* temp = ((HalCbmUnigenEventInterface*) interface)->fEvent;
23 fB = temp->GetB();
24 fPhi = temp->GetPhi();
25 fTotalTracksNo = temp->GetNpa();
26 fTracks->Clear();
27 for (int i = 0; i < fTotalTracksNo; i++) {
28 UParticle* particle = temp->GetParticle(i);
29 TParticlePDG* pdg_part = fPDG->GetParticle(particle->GetPdg());
30 Double_t charge = 0;
31 if (pdg_part) {
32 charge = pdg_part->Charge() / 3.0;
33 }
34 Hal::McTrack* target_track = (Hal::McTrack*) fTracks->ConstructedAt(i);
35 target_track->ResetTrack(i, this);
36 target_track->SetCharge(charge);
37 target_track->SetPdg(particle->GetPdg());
38 if (particle->GetParent() < 0) {
39 target_track->SetPrimary();
40 }
41 else {
42 target_track->SetMotherIndex(particle->GetParent());
43 }
44 target_track->SetMomentum(particle->Px(), particle->Py(), particle->Pz(), particle->E());
45 target_track->SetFreezoutPosition(particle->X(), particle->Y(), particle->Z(), particle->T());
46 target_track->SetStatus(particle->GetStatus());
47 }
48}
49
50void HalCbmUnigenEvent::Clear(Option_t* opt) { Hal::McEvent::Clear(opt); }
51
53
55
56TString HalCbmUnigenEvent::GetFormatName() const { return "UnigenFormat"; }
57
59{
60 Hal::DataManager* manager = Hal::DataManager::Instance();
61 if (manager->CheckBranch("UEvent.")) {
62 return kTRUE;
63 }
64 return kFALSE;
65}
66
67Hal::EventInterface* HalCbmUnigenEvent::CreateInterface() const { return new HalCbmUnigenEventInterface(); }
virtual void Clear(Option_t *opt=" ")
virtual Hal::EventInterface * CreateInterface() const
virtual void Update(Hal::EventInterface *interface)
virtual TString GetFormatName() const
virtual Bool_t ExistInTree() const
Int_t GetNpa() const
Definition UEvent.h:42
Double_t GetPhi() const
Definition UEvent.h:38
UParticle * GetParticle(Int_t index) const
Definition UEvent.cxx:101
Double_t GetB() const
Definition UEvent.h:37
Int_t GetPdg() const
Definition UParticle.h:50
Double_t Py() const
Definition UParticle.h:59
Double_t E() const
Definition UParticle.h:61
Double_t Y() const
Definition UParticle.h:65
Int_t GetParent() const
Definition UParticle.h:52
Double_t Pz() const
Definition UParticle.h:60
Double_t Px() const
Definition UParticle.h:58
Double_t X() const
Definition UParticle.h:64
Double_t Z() const
Definition UParticle.h:66
Int_t GetStatus() const
Definition UParticle.h:51
Double_t T() const
Definition UParticle.h:67