CbmRoot
Loading...
Searching...
No Matches
CbmBeamGenerator.cxx
Go to the documentation of this file.
1/* Copyright (C) 2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese [committer] */
4
11#include "CbmBeamGenerator.h"
12
13#include "CbmEventGenerator.h"
14
15#include "FairMCEventHeader.h"
16#include "FairPrimaryGenerator.h"
17#include "FairRunSim.h"
18#include "Logger.h"
19
20#include "TFile.h"
21#include "TRandom.h"
22#include "TTree.h"
23#include "TVector3.h"
24#include <TDatabasePDG.h>
25#include <TParticlePDG.h>
26
27#include <cassert>
28#include <sstream>
29
30#include "UEvent.h"
31#include "UParticle.h"
32#include "URun.h"
33
34
35// ----- Default constructor --------------------------------------------
37 : FairGenerator("BeamGenerator", "CBM generator")
38 , fP(0.)
39 , fStartZ(0.)
40 , fIon(nullptr)
41{
42}
43// --------------------------------------------------------------------------
44
45
46// ----- Standard constructor --------------------------------------------
47CbmBeamGenerator::CbmBeamGenerator(UInt_t beamZ, UInt_t beamA, UInt_t beamQ, Double_t momentum, Double_t startZ)
48 : FairGenerator("BeamGenerator", "CBM generator")
49 , fP(momentum * Double_t(beamA))
50 , fStartZ(startZ)
51 , fIon(nullptr)
52{
53
54 // --- Create the ion species and add it to the particle list
55 size_t buf_size = 20;
56 char name[buf_size];
57 snprintf(name, buf_size - 1, "Beam_%d_%d_%d", beamZ, beamA, beamQ);
58 fIon = new FairIon(name, beamZ, beamA, beamQ);
59 FairRunSim* run = FairRunSim::Instance();
60 assert(run);
61 run->AddNewIon(fIon);
62
63 // --- Tell the event generator to force the vertex at the given z
64 FairPrimaryGenerator* primGen = run->GetPrimaryGenerator();
65 CbmEventGenerator* evGen = dynamic_cast<CbmEventGenerator*>(primGen);
66 assert(evGen);
67 evGen->ForceVertexAtZ(startZ);
68 evGen->SmearVertexZ(kFALSE);
69}
70// --------------------------------------------------------------------------
71
72
73// ----- Destructor -----------------------------------------------------
75// --------------------------------------------------------------------------
76
77
78// ----- Print info -----------------------------------------------------
79void CbmBeamGenerator::Print(Option_t*) const { LOG(info) << ToString(); }
80// --------------------------------------------------------------------------
81
82
83// ----- Generate input event -------------------------------------------
84Bool_t CbmBeamGenerator::ReadEvent(FairPrimaryGenerator* primGen)
85{
86
87 TParticlePDG* ion = TDatabasePDG::Instance()->GetParticle(fIon->GetName());
88 assert(ion);
89
90 // Note that the beam position in x and y and the beam angle are generated
91 // by CbmEventGenerator. Here, we generate the ion thus in z direction
92 // and at zero coordinates. It will be properly placed and rotated
93 // according to the beam profile specified to CbmEventGenerator.
94 primGen->AddTrack(ion->PdgCode(), 0., 0., fP, 0., 0., 0.);
95
96 return kTRUE;
97}
98// --------------------------------------------------------------------------
99
100
101// ----- Info -----------------------------------------------------------
102std::string CbmBeamGenerator::ToString() const
103{
104
105 std::stringstream ss;
106 ss << GetName() << " ion: " << fIon->GetName() << " Z " << fIon->GetZ() << " A " << fIon->GetA() << " Q "
107 << fIon->GetQ() << " mass " << fIon->GetMass() << ", momentum " << fP << " GeV/u, start Z = " << fStartZ;
108
109 return ss.str();
110}
111// --------------------------------------------------------------------------
112
113
ClassImp(CbmBeamGenerator)
CbmBeamGenerator()
Default constructor (should not be used)
Double_t fStartZ
z coordinate of start point
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
Generate one event (abstract in base class)
FairIon * fIon
Ion type.
std::string ToString() const
Info to string.
virtual ~CbmBeamGenerator()
Destructor.
Double_t fP
Total momentum [GeV].
virtual void Print(Option_t *opt="") const
Print info to logger.
CbmEventGenerator.
void ForceVertexAtZ(Double_t zVertex)
Force event vertex to be at a given z.