Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
CCalMagneticField.cc
이 파일의 문서화 페이지로 가기
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 //
27 // File: CCalMagneticField.cc
28 // Description: User Field class implementation.
30 #include <fstream>
31 
32 #include "CCalMagneticField.hh"
33 #include "CCalutils.hh"
34 #include "G4SystemOfUnits.hh"
35 #include "G4FieldManager.hh"
36 
37 //#define ddebug
38 //#define debug
39 
40 //Constructor and destructor:
41 
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 }
86 
87 
89  if (pos)
90  delete[] pos;
91  if (slope)
92  delete[] slope;
93  if (intercept)
94  delete[] intercept;
95 }
96 
97 
98 // Member functions
99 
100 void CCalMagneticField::MagneticField(const double x[3], double B[3]) const
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 }
133 
134 
136 MagneticField(const CLHEP::Hep3Vector point) const {
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 }
150 
151 
152 void CCalMagneticField::GetFieldValue(const double x[3], double* B) const {
154 }
155 
Float_t x
Definition: compare.C:6
void MagneticField(const double Point[3], double Bfield[3]) const
static const G4double pos
static constexpr double kilogauss
Definition: G4SIunits.hh:271
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
double z() const
virtual void GetFieldValue(const double Point[3], double *Bfield) const
void setX(double)
static constexpr double m
Definition: G4SIunits.hh:129
double G4double
Definition: G4Types.hh:76
CCalMagneticField(const G4String &name)
void setZ(double)
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
std::ostream & tab(std::ostream &)
Definition: CCalutils.cc:89
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
double B(double temperature)
void setY(double)
static constexpr double tesla
Definition: G4SIunits.hh:268