Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagneticFieldModel Class Reference

#include <G4MagneticFieldModel.hh>

Inheritance diagram for G4MagneticFieldModel:

Public Member Functions

 G4MagneticFieldModel (G4int nDataPointsPerHalfExtent=3, Representation representation=Representation::fullArrow, G4int arrow3DLineSegmentsPerCircle=6, const G4VisExtent &extentForField=G4VisExtent(), const std::vector< G4PhysicalVolumesSearchScene::Findings > &pvFindings=std::vector< G4PhysicalVolumesSearchScene::Findings >())
virtual ~G4MagneticFieldModel ()
Public Member Functions inherited from G4VFieldModel
 G4VFieldModel (const G4String &typeOfField, const G4String &symbol="", const G4VisExtent &extentForField=G4VisExtent(), const std::vector< G4PhysicalVolumesSearchScene::Findings > &pvFindings=std::vector< G4PhysicalVolumesSearchScene::Findings >(), G4int nDataPointsPerHalfScene=10, Representation representation=Representation::fullArrow, G4int arrow3DLineSegmentsPerCircle=6)
virtual ~G4VFieldModel ()
virtual void DescribeYourselfTo (G4VGraphicsScene &sceneHandler)
Public Member Functions inherited from G4VModel
 G4VModel (const G4ModelingParameters *=0)
virtual ~G4VModel ()
const G4ModelingParametersGetModelingParameters () const
const G4StringGetType () const
virtual G4String GetCurrentDescription () const
virtual G4String GetCurrentTag () const
const G4VisExtentGetExtent () const
const G4StringGetGlobalDescription () const
const G4StringGetGlobalTag () const
void SetModelingParameters (const G4ModelingParameters *)
void SetExtent (const G4VisExtent &)
void SetType (const G4String &)
void SetGlobalDescription (const G4String &)
void SetGlobalTag (const G4String &)
virtual G4bool Validate (G4bool warn=true)

Protected Member Functions

virtual void GetFieldAtLocation (const G4Field *field, const G4Point3D &position, G4double time, G4Point3D &result) const

Additional Inherited Members

Public Types inherited from G4VFieldModel
enum  Representation { fullArrow , lightArrow }
Static Public Member Functions inherited from G4VModel
static const G4ModelingParametersGetCurrentModelingParameters ()
static void SetCurrentModelingParameters (const G4ModelingParameters *)
Protected Attributes inherited from G4VModel
G4String fType
G4String fGlobalTag
G4String fGlobalDescription
G4VisExtent fExtent
const G4ModelingParametersfpMP

Detailed Description

Definition at line 42 of file G4MagneticFieldModel.hh.

Constructor & Destructor Documentation

◆ G4MagneticFieldModel()

G4MagneticFieldModel::G4MagneticFieldModel ( G4int nDataPointsPerHalfExtent = 3,
Representation representation = Representation::fullArrow,
G4int arrow3DLineSegmentsPerCircle = 6,
const G4VisExtent & extentForField = G4VisExtent(),
const std::vector< G4PhysicalVolumesSearchScene::Findings > & pvFindings = std::vector<G4PhysicalVolumesSearchScene::Findings>() )
inline

Definition at line 46 of file G4MagneticFieldModel.hh.

54 ("Magnetic","B", extentForField, pvFindings,
55 nDataPointsPerHalfExtent, representation, arrow3DLineSegmentsPerCircle)
56 {}
G4VFieldModel(const G4String &typeOfField, const G4String &symbol="", const G4VisExtent &extentForField=G4VisExtent(), const std::vector< G4PhysicalVolumesSearchScene::Findings > &pvFindings=std::vector< G4PhysicalVolumesSearchScene::Findings >(), G4int nDataPointsPerHalfScene=10, Representation representation=Representation::fullArrow, G4int arrow3DLineSegmentsPerCircle=6)

◆ ~G4MagneticFieldModel()

virtual G4MagneticFieldModel::~G4MagneticFieldModel ( )
inlinevirtual

Definition at line 58 of file G4MagneticFieldModel.hh.

58{;}

Member Function Documentation

◆ GetFieldAtLocation()

void G4MagneticFieldModel::GetFieldAtLocation ( const G4Field * field,
const G4Point3D & position,
G4double time,
G4Point3D & result ) const
protectedvirtual

Implements G4VFieldModel.

Definition at line 43 of file G4MagneticFieldModel.cc.

45 {
46 if (!field) return; // No action if no field
47
48 G4double xyzt[4] = { position.x(), position.y(), position.z(), time };
49 G4double BEvals[6] = {0.}; // Field returns {Bx,By,Bz,Ex,Ey,Ez}
50 field->GetFieldValue(xyzt, BEvals);
51
52 result.set(BEvals[0], BEvals[1], BEvals[2]);
53 return;
54}
double G4double
Definition G4Types.hh:83
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
void set(T x1, T y1, T z1)

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