CbmRoot
Loading...
Searching...
No Matches
CbmEventGenerator Class Reference

CbmEventGenerator. More...

#include <CbmEventGenerator.h>

Inheritance diagram for CbmEventGenerator:
[legend]
Collaboration diagram for CbmEventGenerator:
[legend]

Public Member Functions

 CbmEventGenerator ()
 Default constructor.
 
virtual ~CbmEventGenerator ()
 Destructor.
 
virtual void AddTrack (Int_t pdgid, Double_t px, Double_t py, Double_t pz, Double_t vx, Double_t vy, Double_t vz, Int_t parent=-1, Bool_t wanttracking=true, Double_t e=-9e9, Double_t tof=0., Double_t weight=0., TMCProcess proc=kPPrimary)
 Add a track to the stack to be transported.
 
void ForceVertexAtZ (Double_t zVertex)
 Force event vertex to be at a given z.
 
void ForceVertexInTarget (Bool_t choice=kTRUE)
 Enable or disable forcing the vertex to be in the target.
 
const CbmBeamProfileGetBeamProfile ()
 Beam profile.
 
virtual Bool_t GenerateEvent (FairGenericStack *stack)
 Generate the input event.
 
virtual void Print (Option_t *opt="") const
 Log output.
 
void SetBeamAngle (Double_t meanThetaX, Double_t meanThetaY, Double_t sigmaThetaX=-1., Double_t sigmaThetaY=-1.)
 Set the beam angle in the focal plane.
 
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.
 
void SetTarget (std::shared_ptr< const CbmTarget > target)
 Set target properties.
 

Private Member Functions

virtual void MakeBeamAngle ()
 Generate beam angle.
 
virtual void MakeVertex ()
 Generate event vertex position.
 
void MakeVertexAtZ ()
 Generate event vertex position at a given z.
 
virtual void MakeVertexInFocalPlane ()
 Generate event vertex position in the beam focal plane.
 
virtual void MakeVertexInTarget ()
 Generate event vertex position in the target.
 
 ClassDef (CbmEventGenerator, 2)
 

Private Attributes

CbmBeamProfile fBeamProfile
 Beam properties.
 
std::shared_ptr< const CbmTargetfTarget
 
Bool_t fForceVertexInTarget
 Target properties.
 
Bool_t fForceVertexAtZ
 If set, vertex must be at given z.
 
Double_t fVertexZ
 forced z coordinate of event vertex
 

Detailed Description

CbmEventGenerator.

Author
Volker Friese v.fri.nosp@m.ese@.nosp@m.gsi.d.nosp@m.e
Since
29 July 2019
Date
29 July 2019

The EventGenerator defines the primary particles as input to the transport simulation. Several event generator objects can be registered to the EventGenerator, each generating primary particles. The EventGenerator generates the common event vertex and rotation according to the user specification.

CbmEventGenerator derives from FairPrimaryGenerator. It re-implements the methods MakeBeam and MakeVertex such as to ensure that the generated vertex falls into the target volume. It also re-implements te method AddTrack to properly implement the rotation following the beam direction.

Definition at line 41 of file CbmEventGenerator.h.

Constructor & Destructor Documentation

◆ CbmEventGenerator()

CbmEventGenerator::CbmEventGenerator ( )

Default constructor.

Definition at line 31 of file CbmEventGenerator.cxx.

References fBeamProfile, fForceVertexAtZ, fForceVertexInTarget, fTarget, and fVertexZ.

Referenced by ClassDef().

◆ ~CbmEventGenerator()

CbmEventGenerator::~CbmEventGenerator ( )
virtual

Destructor.

Definition at line 47 of file CbmEventGenerator.cxx.

Member Function Documentation

◆ AddTrack()

void CbmEventGenerator::AddTrack ( Int_t pdgid,
Double_t px,
Double_t py,
Double_t pz,
Double_t vx,
Double_t vy,
Double_t vz,
Int_t parent = -1,
Bool_t wanttracking = true,
Double_t e = -9e9,
Double_t tof = 0.,
Double_t weight = 0.,
TMCProcess proc = kPPrimary )
virtual

Add a track to the stack to be transported.

Definition at line 54 of file CbmEventGenerator.cxx.

◆ ClassDef()

CbmEventGenerator::ClassDef ( CbmEventGenerator ,
2  )
private

References CbmEventGenerator().

◆ ForceVertexAtZ()

void CbmEventGenerator::ForceVertexAtZ ( Double_t zVertex)

Force event vertex to be at a given z.

Parameters
zVertexz coordinate of event vertex

The event vertex will be sampled in x and y from the specified beam properties (profile in focal plane and angular distribution),

Definition at line 144 of file CbmEventGenerator.cxx.

References fForceVertexAtZ, fForceVertexInTarget, and fVertexZ.

Referenced by CbmBeamGenerator::CbmBeamGenerator().

◆ ForceVertexInTarget()

void CbmEventGenerator::ForceVertexInTarget ( Bool_t choice = kTRUE)
inline

Enable or disable forcing the vertex to be in the target.

Parameters
choiceIf true, the vertex will be generated in the target

Definition at line 70 of file CbmEventGenerator.h.

References fForceVertexInTarget.

◆ GenerateEvent()

Bool_t CbmEventGenerator::GenerateEvent ( FairGenericStack * stack)
virtual

Generate the input event.

Parameters
stackPointer to stack object

This is the main functionality, invoked from FairRunSim. It fills the stack at the beginning of each event with the primary particles provided by the generator instances. All primary track momenta are rotated according to the beam direction (in x-z and y-z) and and by the event plane (in x-y).

Re-implemented from base-class FairPrimaryGenerator.

Definition at line 154 of file CbmEventGenerator.cxx.

References MakeVertex().

◆ GetBeamProfile()

const CbmBeamProfile & CbmEventGenerator::GetBeamProfile ( )
inline

Beam profile.

Returns
Reference to beam profile object

Definition at line 76 of file CbmEventGenerator.h.

References fBeamProfile.

◆ MakeBeamAngle()

virtual void CbmEventGenerator::MakeBeamAngle ( )
inlineprivatevirtual

Generate beam angle.

          Will be called from FairPrimaryGenerator::GenerateEvent().
          The method is re-implemented here as empty. The beam angle
          will be generated along with the beam position in the method
          MakeVertex().

Definition at line 157 of file CbmEventGenerator.h.

◆ MakeVertex()

void CbmEventGenerator::MakeVertex ( )
privatevirtual

Generate event vertex position.

          Will be called from FairPrimaryGenerator::GenerateEvent().
          Beam position and angles in the focal plane are sampled
          from the specified distributions. There are three options to generate
          the event vertex:
          1. If a vertex z position is specified by a call to ForceVertexAtZ,

the event vertex is the beam extrapolation to the specified z.

  1. Else, if a target is specified, the event vertex is sampled from the beam trajectory within the target - flat distributions between entry and exit point, or always in the target centre plane, depending on the choice set with SetSmearVertexZ().
  2. Else, the event vertex is the beam position in the focal plane.

Definition at line 216 of file CbmEventGenerator.cxx.

References fForceVertexAtZ, fForceVertexInTarget, fTarget, MakeVertexAtZ(), MakeVertexInFocalPlane(), and MakeVertexInTarget().

Referenced by GenerateEvent().

◆ MakeVertexAtZ()

void CbmEventGenerator::MakeVertexAtZ ( )
private

Generate event vertex position at a given z.

Will be used if ForceVertexAtZ was called.

A point in the plane z = zVertex

Normal on the plane z = const.

Definition at line 232 of file CbmEventGenerator.cxx.

References fBeamProfile, and fVertexZ.

Referenced by MakeVertex().

◆ MakeVertexInFocalPlane()

void CbmEventGenerator::MakeVertexInFocalPlane ( )
privatevirtual

Generate event vertex position in the beam focal plane.

          Will be used if no target was specified.

Definition at line 248 of file CbmEventGenerator.cxx.

References fBeamProfile.

Referenced by MakeVertex().

◆ MakeVertexInTarget()

void CbmEventGenerator::MakeVertexInTarget ( )
privatevirtual

Generate event vertex position in the target.

          Will be called if a target was specified.

Definition at line 262 of file CbmEventGenerator.cxx.

References fBeamProfile, and fTarget.

Referenced by MakeVertex().

◆ Print()

void CbmEventGenerator::Print ( Option_t * opt = "") const
virtual

Log output.

Definition at line 319 of file CbmEventGenerator.cxx.

References fBeamProfile, and fTarget.

◆ SetBeamAngle()

void CbmEventGenerator::SetBeamAngle ( Double_t meanThetaX,
Double_t meanThetaY,
Double_t sigmaThetaX = -1.,
Double_t sigmaThetaY = -1. )

Set the beam angle in the focal plane.

Parameters
meanThetaXMean angle in x-z [rad]
meanThetaYMean angle in y-z [rad]
sigmaThetaXRMS of beam angle in x-z [rad]
sigmaThetaYRMS of beam angle in y-z [rad]

The beam angles will be sampled from Gaussian distributions with the specified mean values and RMS. If the latter are negative, the beam angles will be fixed to their mean values.

Default is (0., 0.), no sampling.

Note: Re-implements the non-virtual method in FairPrimaryGenerator.

Definition at line 339 of file CbmEventGenerator.cxx.

References fBeamProfile.

◆ SetBeamPosition()

void CbmEventGenerator::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.

Parameters
meanXMean x position [cm]
meanYMean y position [cm|
sigmaXRMS of x position [cm]
sigmaYRMS of y position [cm]
zFz position of focal plane [cm]

If the beam widths in x and/or y are positive, the beam position will be sampled randomly for each event from a Gauss distribution. Otherwise, it is fixed for all events.

Default is (0.,0.,0.) for the position, no sampling.

Definition at line 348 of file CbmEventGenerator.cxx.

References fBeamProfile.

◆ SetTarget()

void CbmEventGenerator::SetTarget ( std::shared_ptr< const CbmTarget > target)
inline

Set target properties.

Parameters
targetPointer to CbmTarget instance

When a target is specified, the event vertex will be sampled by a constant distribution along the (straight) trajectory of the beam inside the target. If SmearVertexZ(kFALSE) is set afterwards, the vertex will always be in the target centre plane.

Definition at line 139 of file CbmEventGenerator.h.

References fTarget.

Member Data Documentation

◆ fBeamProfile

CbmBeamProfile CbmEventGenerator::fBeamProfile
private

◆ fForceVertexAtZ

Bool_t CbmEventGenerator::fForceVertexAtZ
private

If set, vertex must be at given z.

Definition at line 146 of file CbmEventGenerator.h.

Referenced by CbmEventGenerator(), ForceVertexAtZ(), and MakeVertex().

◆ fForceVertexInTarget

Bool_t CbmEventGenerator::fForceVertexInTarget
private

Target properties.

If set, vertex must be in target

Definition at line 145 of file CbmEventGenerator.h.

Referenced by CbmEventGenerator(), ForceVertexAtZ(), ForceVertexInTarget(), and MakeVertex().

◆ fTarget

std::shared_ptr<const CbmTarget> CbmEventGenerator::fTarget
private

◆ fVertexZ

Double_t CbmEventGenerator::fVertexZ
private

forced z coordinate of event vertex

Definition at line 147 of file CbmEventGenerator.h.

Referenced by CbmEventGenerator(), ForceVertexAtZ(), and MakeVertexAtZ().


The documentation for this class was generated from the following files: