CbmRoot
|
#include <CbmEventGenerator.h>
Public Member Functions | |
CbmEventGenerator () | |
Default constructor More... | |
virtual | ~CbmEventGenerator () |
Destructor More... | |
void | ForceVertexAtZ (Double_t zVertex) |
Force event vertex to be at a given z. More... | |
void | ForceVertexInTarget (Bool_t choice=kTRUE) |
Enable or disable forcing the vertex to be in the target. More... | |
const CbmBeamProfile & | GetBeamProfile () |
Beam profile. More... | |
virtual Bool_t | GenerateEvent (FairGenericStack *stack) |
Generate the input event. More... | |
virtual void | Print (Option_t *opt="") const |
Log output. More... | |
void | SetBeamAngle (Double_t meanThetaX, Double_t meanThetaY, Double_t sigmaThetaX=-1., Double_t sigmaThetaY=-1.) |
Set the beam angle in the focal plane. More... | |
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. More... | |
void | SetTarget (std::shared_ptr< const CbmTarget > target) |
Set target properties. More... | |
Private Member Functions | |
virtual void | MakeBeamAngle () |
Generate beam angle. More... | |
virtual void | MakeVertex () |
Generate event vertex position. More... | |
void | MakeVertexAtZ () |
Generate event vertex position at a given z. More... | |
virtual void | MakeVertexInFocalPlane () |
Generate event vertex position in the beam focal plane. More... | |
virtual void | MakeVertexInTarget () |
Generate event vertex position in the target. More... | |
ClassDef (CbmEventGenerator, 2) | |
Private Attributes | |
CbmBeamProfile | fBeamProfile |
Beam properties. More... | |
std::shared_ptr< const CbmTarget > | fTarget |
Bool_t | fForceVertexInTarget |
Target properties. More... | |
Bool_t | fForceVertexAtZ |
If set, vertex must be at given z. More... | |
Double_t | fVertexZ |
forced z coordinate of event vertex More... | |
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.
Definition at line 34 of file CbmEventGenerator.h.
CbmEventGenerator::CbmEventGenerator | ( | ) |
Default constructor
Definition at line 19 of file CbmEventGenerator.cxx.
|
virtual |
Destructor
Definition at line 34 of file CbmEventGenerator.cxx.
|
private |
void CbmEventGenerator::ForceVertexAtZ | ( | Double_t | zVertex | ) |
Force event vertex to be at a given z.
zVertex | z 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 39 of file CbmEventGenerator.cxx.
References fForceVertexAtZ, fForceVertexInTarget, and fVertexZ.
Referenced by CbmBeamGenerator::CbmBeamGenerator(), and CbmTransport::ForceVertexAtZ().
|
inline |
Enable or disable forcing the vertex to be in the target.
choice | If true, the vertex will be generated in the target |
Definition at line 57 of file CbmEventGenerator.h.
References fForceVertexInTarget.
Referenced by CbmTransport::ForceVertexInTarget().
|
virtual |
Generate the input event.
stack | Pointer 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 48 of file CbmEventGenerator.cxx.
References MakeVertex().
|
inline |
Beam profile.
Definition at line 65 of file CbmEventGenerator.h.
References fBeamProfile.
Referenced by CbmTransport::InitEventGenerator().
|
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 153 of file CbmEventGenerator.h.
|
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.
Definition at line 113 of file CbmEventGenerator.cxx.
References fForceVertexAtZ, fForceVertexInTarget, fTarget, MakeVertexAtZ(), MakeVertexInFocalPlane(), and MakeVertexInTarget().
Referenced by GenerateEvent().
|
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 126 of file CbmEventGenerator.cxx.
References fBeamProfile, fVertexZ, and CbmBeamProfile::GenerateBeam().
Referenced by MakeVertex().
|
privatevirtual |
Generate event vertex position in the beam focal plane.
Will be used if no target was specified.
Definition at line 141 of file CbmEventGenerator.cxx.
References fBeamProfile, and CbmBeamProfile::GenerateBeam().
Referenced by MakeVertex().
|
privatevirtual |
Generate event vertex position in the target.
Will be called if a target was specified.
Definition at line 154 of file CbmEventGenerator.cxx.
References fBeamProfile, fTarget, and CbmBeamProfile::GenerateBeam().
Referenced by MakeVertex().
|
virtual |
Log output.
Definition at line 212 of file CbmEventGenerator.cxx.
References fBeamProfile, fTarget, and CbmBeamProfile::ToString().
Referenced by CbmTransport::InitEventGenerator().
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.
meanThetaX | Mean angle in x-z [rad] |
meanThetaY | Mean angle in y-z [rad] |
sigmaThetaX | RMS of beam angle in x-z [rad] |
sigmaThetaY | RMS 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 231 of file CbmEventGenerator.cxx.
References fBeamProfile, and CbmBeamProfile::SetAngle().
Referenced by CbmTransport::SetBeamAngle().
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.
meanX | Mean x position [cm] |
meanY | Mean y position [cm| |
sigmaX | RMS of x position [cm] |
sigmaY | RMS of y position [cm] |
zF | z 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 241 of file CbmEventGenerator.cxx.
References fBeamProfile, and CbmBeamProfile::SetPosition().
Referenced by CbmTransport::SetBeamPosition().
|
inline |
Set target properties.
target | Pointer 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 135 of file CbmEventGenerator.h.
References fTarget.
Referenced by CbmTransport::InitEventGenerator().
|
private |
Beam properties.
Definition at line 139 of file CbmEventGenerator.h.
Referenced by GetBeamProfile(), MakeVertexAtZ(), MakeVertexInFocalPlane(), MakeVertexInTarget(), Print(), SetBeamAngle(), and SetBeamPosition().
|
private |
If set, vertex must be at given z.
Definition at line 142 of file CbmEventGenerator.h.
Referenced by ForceVertexAtZ(), and MakeVertex().
|
private |
Target properties.
If set, vertex must be in target
Definition at line 141 of file CbmEventGenerator.h.
Referenced by ForceVertexAtZ(), ForceVertexInTarget(), and MakeVertex().
|
private |
Definition at line 140 of file CbmEventGenerator.h.
Referenced by MakeVertex(), MakeVertexInTarget(), Print(), and SetTarget().
|
private |
forced z coordinate of event vertex
Definition at line 143 of file CbmEventGenerator.h.
Referenced by ForceVertexAtZ(), and MakeVertexAtZ().