Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CCalMagneticField Class Reference

#include <CCalMagneticField.hh>

Inheritance diagram for CCalMagneticField:
Collaboration diagram for CCalMagneticField:

Public Member Functions

 CCalMagneticField (const G4String &name)
 
 ~CCalMagneticField ()
 
void MagneticField (const double Point[3], double Bfield[3]) const
 
CLHEP::Hep3Vector MagneticField (const CLHEP::Hep3Vector Point) const
 
virtual void GetFieldValue (const double Point[3], double *Bfield) const
 
G4double GetConstantFieldvalue () const
 
- Public Member Functions inherited from G4MagneticField
 G4MagneticField ()
 
virtual ~G4MagneticField ()
 
 G4MagneticField (const G4MagneticField &r)
 
G4MagneticFieldoperator= (const G4MagneticField &p)
 
G4bool DoesFieldChangeEnergy () const
 
virtual void GetFieldValue (const G4double Point[4], G4double *Bfield) const =0
 
- Public Member Functions inherited from G4ElectroMagneticField
 G4ElectroMagneticField ()
 
virtual ~G4ElectroMagneticField ()
 
 G4ElectroMagneticField (const G4ElectroMagneticField &r)
 
G4ElectroMagneticFieldoperator= (const G4ElectroMagneticField &p)
 
- Public Member Functions inherited from G4Field
 G4Field (G4bool gravityOn=false)
 
 G4Field (const G4Field &)
 
virtual ~G4Field ()
 
G4Fieldoperator= (const G4Field &p)
 
G4bool IsGravityActive () const
 
void SetGravityActive (G4bool OnOffFlag)
 
virtual G4FieldClone () const
 

Protected Member Functions

G4FieldManagerGetGlobalFieldManager ()
 

Detailed Description

Definition at line 38 of file CCalMagneticField.hh.

Constructor & Destructor Documentation

CCalMagneticField::CCalMagneticField ( const G4String name)

Definition at line 42 of file CCalMagneticField.cc.

42  :
43  fval(0), pos(0), slope(0), intercept(0)
44 {
45 #ifdef debug
46  fVerbosity = 1;
47 #else
48  fVerbosity = 0;
49 #endif
50 
51  //Let's open the file
52  G4cout << " ==> Opening file " << filename << " to read magnetic field..."
53  << G4endl;
54  G4String pathName = getenv("CCAL_GLOBALPATH");
55  std::ifstream is;
56  bool ok = openGeomFile(is, pathName, filename);
57 
58  if (ok) {
59  findDO(is, G4String("FLDM"));
60  is >> fval >> npts >> xoff;
61 
62  if (fVerbosity)
63  G4cout << "Field value " << fval << " # points " << npts <<
64  " offset in x "
65  << xoff*mm << G4endl;
66 
67  if (npts > 0) {
68  pos = new G4double[npts];
69  slope = new G4double[npts];
70  intercept = new G4double[npts];
71 
72  for (G4int i = 0; i < npts; i++) {
73  is >> pos[i] >> slope[i] >> intercept[i];
74  if (fVerbosity)
75  G4cout << tab << "Position " << i << " " << pos[i] << " Slope "
76  << slope[i] << " Intercept " << intercept[i] << G4endl;
77  }
78  }
79 
81  // Close the file
82  G4cout << " ==> Closing file " << filename << G4endl;
83  is.close();
84  }
85 }
static constexpr double mm
Definition: G4SIunits.hh:115
bool openGeomFile(std::ifstream &is, const G4String &pathname, const G4String &filename)
Definition: CCalutils.cc:116
std::ifstream & findDO(std::ifstream &, const G4String &)
Definition: CCalutils.cc:72
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
std::ostream & tab(std::ostream &)
Definition: CCalutils.cc:89
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos

Here is the call graph for this function:

CCalMagneticField::~CCalMagneticField ( )

Definition at line 88 of file CCalMagneticField.cc.

88  {
89  if (pos)
90  delete[] pos;
91  if (slope)
92  delete[] slope;
93  if (intercept)
94  delete[] intercept;
95 }
static const G4double pos

Member Function Documentation

G4double CCalMagneticField::GetConstantFieldvalue ( ) const
inline

Definition at line 47 of file CCalMagneticField.hh.

47 {return fval;}

Here is the caller graph for this function:

void CCalMagneticField::GetFieldValue ( const double  Point[3],
double *  Bfield 
) const
virtual

Definition at line 152 of file CCalMagneticField.cc.

152  {
154 }
tuple x
Definition: test.py:50
double B(double temperature)
void MagneticField(const double Point[3], double Bfield[3]) const

Here is the call graph for this function:

G4FieldManager* CCalMagneticField::GetGlobalFieldManager ( )
protected
void CCalMagneticField::MagneticField ( const double  Point[3],
double  Bfield[3] 
) const

Definition at line 100 of file CCalMagneticField.cc.

101 {
102  G4int i=0;
103  for (i=0; i<2; i++) {
104  B[i] = 0*kilogauss;
105  }
106 
107  G4double m1=0;
108  G4double c1=0;
109  G4double xnew = x[0]/mm + xoff;
110  if (npts > 0) {
111  for (i=0; i<npts; i++) {
112  if (xnew > pos[i]*mm) {
113  m1 = slope[i];
114  c1 = intercept[i];
115  }
116  }
117  }
118  G4double scor = c1 + m*xnew;
119  if (scor < 0.) scor = 0.;
120  if (scor > 1.) scor = 1.0;
121 
122  B[2] = scor*fval*kilogauss;
123  if (fVerbosity)
124  {
125 
126  G4cout << "Field at x: " << x[0]/mm << "mm (" << xnew << ") = " <<
127  B[2]/tesla
128  << "T (m = " << m1 << ", c = " <<
129  c1 << ", scale = " << scor << ")"
130  << G4endl;
131  }
132 }
static constexpr double tesla
Definition: G4SIunits.hh:268
static constexpr double mm
Definition: G4SIunits.hh:115
tuple x
Definition: test.py:50
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double kilogauss
Definition: G4SIunits.hh:271
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos
tuple c1
Definition: plottest35.py:14

Here is the caller graph for this function:

CLHEP::Hep3Vector CCalMagneticField::MagneticField ( const CLHEP::Hep3Vector  Point) const

Definition at line 136 of file CCalMagneticField.cc.

136  {
137 
138  G4double x[3],B[3];
140 
141  x[0] = point.x();
142  x[1] = point.y();
143  x[2] = point.z();
145  v.setX(B[0]);
146  v.setY(B[1]);
147  v.setZ(B[2]);
148  return v;
149 }
tuple x
Definition: test.py:50
double B(double temperature)
void setY(double)
void setZ(double)
void setX(double)
tuple v
Definition: test.py:18
void MagneticField(const double Point[3], double Bfield[3]) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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