CbmRoot
Loading...
Searching...
No Matches
CbmEventGenerator.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese [committer] */
4
10#include "CbmEventGenerator.h"
11
12#include "CbmBeam.h"
13
14#include <FairGenericStack.h>
15#include <FairMCEventHeader.h>
16#include <Logger.h>
17
18#include <TMath.h>
19#include <TRandom.h>
20#include <TVector3.h>
21
22#include <cassert>
23
24// ----- Constructor ---------------------------------------------------
26 : FairPrimaryGenerator()
27 , fBeamProfile()
28 , fTarget()
29 , fForceVertexInTarget(kTRUE)
30 , fForceVertexAtZ(kFALSE)
31 , fVertexZ(0.)
32{
33 // Mother class members
34 fName = "EventGenerator";
35 fBeamAngle = kTRUE;
36}
37// -------------------------------------------------------------------------
38
39
40// ----- Destructor ----------------------------------------------------
42// -------------------------------------------------------------------------
43
44
45// ----- Force vertex at a given z position ----------------------------
47{
48 fForceVertexAtZ = kTRUE;
49 fForceVertexInTarget = kFALSE;
50 fVertexZ = zVertex;
51}
52// -------------------------------------------------------------------------
53
54
55// ----- Generate the input event --------------------------------------
56Bool_t CbmEventGenerator::GenerateEvent(FairGenericStack* stack)
57{
58
59 std::cout << std::endl;
60 LOG(info) << GetName() << ": Generate event " << ++fEventNr;
61
62 // An MCEventHeader must be present
63 assert(fEvent);
64
65 // Set the stack
66 assert(stack);
67 fStack = stack;
68
69 // Reset MCEventHeader and track counter
70 fEvent->Reset();
71 fNTracks = 0;
72
73 // Generate the event vertex and the beam angle
74 MakeVertex();
75
76 // Generate the event plane angle
77 if (fEventPlane) {
78 MakeEventPlane();
79 LOG(info) << GetName() << ": Rotate event by " << fPhi << " rad";
80 }
81
82 // Call the registered generator classes
83 fListIter->Reset();
84 TObject* item = nullptr;
85 while ((item = fListIter->Next())) {
86 FairGenerator* generator = dynamic_cast<FairGenerator*>(item);
87 assert(generator);
88 fMCIndexOffset = fNTracks; // Number of tracks before call of generator
89 if (!generator->ReadEvent(this)) {
90 LOG(info) << GetName() << ": ReadEvent failed for generator " << generator->GetName();
91 return kFALSE;
92 }
93 }
94
95 // Update event header
96 if (!fEvent->IsSet()) fEvent->SetEventID(fEventNr);
97 fEvent->SetNPrim(fNTracks);
98 fEvent->SetVertex(fVertex);
99 fEvent->SetRotX(fBeamAngleX);
100 fEvent->SetRotY(fBeamAngleY);
101 Double_t phiGen = fEvent->GetRotZ(); // from concrete generators
102 fEvent->SetRotZ(phiGen + fPhi); // total event plane angle
103 // Event info to screen
104 LOG(info) << GetName() << ": Event ID " << fEvent->GetEventID() << ", tracks " << fNTracks << ", vertex ("
105 << fVertex.X() << ", " << fVertex.Y() << ", " << fVertex.Z() << ") cm";
106 LOG(info) << GetName() << ": Beam angle (" << fBeamAngleX << ", " << fBeamAngleY << ") rad, event plane angle "
107 << fPhi << " rad (total " << fEvent->GetRotZ() << " rad) change";
108
109 // Update global track counter (who knows what it's good for ...)
110 fTotPrim += fNTracks;
111
112 return kTRUE;
113}
114// -------------------------------------------------------------------------
115
116
117// ----- Generate event vertex -----------------------------------------
127// -------------------------------------------------------------------------
128
129
130// ----- Generate event vertex in the beam focal plane -----------------
132{
133
134 std::unique_ptr<CbmBeam> beam;
135 beam = fBeamProfile.GenerateBeam();
136 TVector3 point(0., 0., fVertexZ);
137 TVector3 norm(0, 0., 1.);
138 fVertex = beam->ExtrapolateToPlane(point, norm);
139 fBeamAngleX = beam->GetThetaX();
140 fBeamAngleY = beam->GetThetaY();
141 fBeamDirection = beam->GetDirection();
142}
143// -------------------------------------------------------------------------
144
145
146// ----- Generate event vertex in the beam focal plane -----------------
148{
149
150 std::unique_ptr<CbmBeam> beam;
151 beam = fBeamProfile.GenerateBeam();
152 fVertex = beam->GetPosition();
153 fBeamAngleX = beam->GetThetaX();
154 fBeamAngleY = beam->GetThetaY();
155 fBeamDirection = beam->GetDirection();
156}
157// -------------------------------------------------------------------------
158
159
160// ----- Generate event vertex -----------------------------------------
162{
163
164 assert(fTarget);
165
166 // Target surface centres and normal vector
167 TVector3 surf1 = fTarget->GetSurfaceCentreUp();
168 TVector3 surf2 = fTarget->GetSurfaceCentreDown();
169 TVector3 norm = fTarget->GetNormal();
170
171 // Declare some vectors
172 TVector3 point1; // Intersection of beam with first target surface
173 TVector3 point2; // Intersection of beam with second target surface
174
175 // It is possible that a generated beam trajectory misses the target.
176 // In that case, the beam will be re-sampled until it hits the target.
177 Bool_t isInTarget = kFALSE;
178 Int_t nSamples = 0;
179 std::unique_ptr<CbmBeam> beam;
180 while (!isInTarget) {
181
182 // Abort if too many beam samples
183 if (nSamples > 1000.)
184 LOG(fatal) << GetName() << ": Aborting after " << nSamples << " beam samplings. Adjust beam and target.";
185
186 // Sample a beam trajectory
187 beam = fBeamProfile.GenerateBeam();
188 nSamples++;
189
190 // Intersections of beam trajectory with target surfaces
191 point1 = beam->ExtrapolateToPlane(surf1, norm);
192 point2 = beam->ExtrapolateToPlane(surf2, norm);
193
194 // Check whether both intersections are inside the target
195 // (i.e., the beam crosses both target surfaces)
196 Bool_t check1 = ((point1 - surf1).Mag() < 0.5 * fTarget->GetDiameter());
197 Bool_t check2 = ((point2 - surf2).Mag() < 0.5 * fTarget->GetDiameter());
198 isInTarget = check1 && check2;
199
200 } //? Beam not in target
201
202 LOG(debug) << beam->ToString() << ", generated after " << nSamples << (nSamples > 1 ? " samplings: " : " sampling");
203
204 // Sample a point between entry and exit of beam in target
205 Double_t scale = 0.5;
206 if (fSmearVertexZ) scale = gRandom->Uniform();
207 fVertex = point1 + scale * (point2 - point1);
208
209 // Set the beam angles (member variables of base class)
210 fBeamAngleX = beam->GetThetaX();
211 fBeamAngleY = beam->GetThetaY();
212 fBeamDirection = beam->GetDirection();
213}
214// -------------------------------------------------------------------------
215
216
217// ----- Info ----------------------------------------------------------
218void CbmEventGenerator::Print(Option_t*) const
219{
220
221 LOG(info) << fBeamProfile.ToString();
222 if (fTarget) LOG(info) << fTarget->ToString();
223 LOG(info) << "Vertex smearing along beam " << (fSmearVertexZ ? "ON" : "OFF");
224 if (fEventPlane) LOG(info) << "Random event plane angle between " << fPhiMin << " and " << fPhiMax << " rad";
225 else
226 LOG(info) << "Fixed event plane angle = 0";
227 LOG(info) << "Number of generators " << fGenList->GetEntries();
228 for (Int_t iGen = 0; iGen < fGenList->GetEntries(); iGen++) {
229 fGenList->At(iGen)->Print();
230 }
231}
232// -------------------------------------------------------------------------
233
234
235// ----- Set beam angle ------------------------------------------------
236void CbmEventGenerator::SetBeamAngle(Double_t meanThetaX, Double_t meanThetaY, Double_t sigmaThetaX,
237 Double_t sigmaThetaY)
238{
239 fBeamProfile.SetAngle(meanThetaX, meanThetaY, sigmaThetaX, sigmaThetaY);
240}
241// -------------------------------------------------------------------------
242
243
244// ----- Set beam position ---------------------------------------------
245void CbmEventGenerator::SetBeamPosition(Double_t meanX, Double_t meanY, Double_t sigmaX, Double_t sigmaY, Double_t zF)
246{
247 fBeamProfile.SetPosition(meanX, meanY, sigmaX, sigmaY, zF);
248}
249// -------------------------------------------------------------------------
250
251
ClassImp(CbmConverterManager)
void SetAngle(Double_t x0, Double_t y0, Double_t sigmaX=-1., Double_t sigmaY=-1.)
Set the parameters for the beam angle distribution.
void SetPosition(Double_t x0, Double_t y0, Double_t sigmaX=-1., Double_t sigmaY=-1., Double_t zF=0.)
Set the parameters for the beam position distribution.
std::string ToString() const
Info to string.
std::unique_ptr< CbmBeam > GenerateBeam()
Generate a beam trajectory.
CbmEventGenerator.
virtual void Print(Option_t *opt="") const
Log output.
virtual Bool_t GenerateEvent(FairGenericStack *stack)
Generate the input event.
Double_t fVertexZ
forced z coordinate of event vertex
void SetBeamAngle(Double_t meanThetaX, Double_t meanThetaY, Double_t sigmaThetaX=-1., Double_t sigmaThetaY=-1.)
Set the beam angle in the focal plane.
CbmEventGenerator()
Default constructor
void MakeVertexAtZ()
Generate event vertex position at a given z.
Bool_t fForceVertexInTarget
Target properties.
virtual ~CbmEventGenerator()
Destructor
Bool_t fForceVertexAtZ
If set, vertex must be at given z.
void ForceVertexAtZ(Double_t zVertex)
Force event vertex to be at a given z.
virtual void MakeVertexInFocalPlane()
Generate event vertex position in the beam focal plane.
virtual void MakeVertex()
Generate event vertex position.
CbmBeamProfile fBeamProfile
Beam properties.
void SetBeamPosition(Double_t meanX, Double_t meanY, Double_t sigmaX=-1., Double_t sigmaY=-1., Double_t zF=0.)
Set the beam position in the focal plane.
std::shared_ptr< const CbmTarget > fTarget
virtual void MakeVertexInTarget()
Generate event vertex position in the target.