Geant4-11
G4TUniformMagneticField.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27#ifndef G4UniformMagneticField_HH
28#define G4UniformMagneticField_HH
29
30#include "G4Types.hh"
31#include "G4ThreeVector.hh"
32#include "G4MagneticField.hh"
33
35{
36 public: // with description
37
39 // A field with value equal to FieldVector.
40 {
41 fFieldComponents[0] = FieldVector.x();
42 fFieldComponents[1] = FieldVector.y();
43 fFieldComponents[2] = FieldVector.z();
44 }
45
46
48 G4double vTheta,
49 G4double vPhi )
50 {
51 if ( (vField<0) || (vTheta<0) || (vTheta>pi) || (vPhi<0) || (vPhi>twopi) )
52 {
53 G4Exception("G4TUniformMagneticField::G4TUniformMagneticField()",
54 "GeomField0002", FatalException, "Invalid parameters.") ;
55 }
56 fFieldComponents[0] = vField*std::sin(vTheta)*std::cos(vPhi) ;
57 fFieldComponents[1] = vField*std::sin(vTheta)*std::sin(vPhi) ;
58 fFieldComponents[2] = vField*std::cos(vTheta) ;
59 }
60
62
65 {
66 for (G4int i=0; i<3; ++i)
68 }
69
71 // Copy constructor and assignment operator.
72 {
73 if (&p == this) return *this;
74 for (G4int i=0; i<3; ++i)
76 return *this;
77 }
78
79 inline void GetFieldValue(const G4double yTrack[4],
80 G4double *B) const
81 {
82 B[0]= fFieldComponents[0] ;
83 B[1]= fFieldComponents[1] ;
84 B[2]= fFieldComponents[2] ;
85 }
86
87 void SetFieldValue(const G4ThreeVector& newFieldVector)
88 {
89 fFieldComponents[0] = newFieldVector.x();
90 fFieldComponents[1] = newFieldVector.y();
91 fFieldComponents[2] = newFieldVector.z();
92 }
93
95 {
99 return B;
100 }
101 // Return the field value
102
104 {
106 this->fFieldComponents[1],
107 this->fFieldComponents[2]) );
108 }
109
110 private:
112};
113
114#endif
G4double B(G4double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double pi
Definition: G4SIunits.hh:55
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
G4TUniformMagneticField(const G4ThreeVector &FieldVector)
G4TUniformMagneticField & operator=(const G4TUniformMagneticField &p)
G4TUniformMagneticField(G4double vField, G4double vTheta, G4double vPhi)
void GetFieldValue(const G4double yTrack[4], G4double *B) const
G4TUniformMagneticField(const G4TUniformMagneticField &p)
void SetFieldValue(const G4ThreeVector &newFieldVector)
virtual G4TUniformMagneticField * Clone() const
G4ThreeVector GetConstantFieldValue() const