Geant4-11
Public Member Functions | Private Attributes
G4ChannelingECHARM Class Reference

#include <G4ChannelingECHARM.hh>

Public Member Functions

 G4ChannelingECHARM (const G4String &, G4double)
 
G4double GetEC (G4ThreeVector &)
 
G4double GetIntSp (G4int index)
 
G4double GetMax ()
 
G4double GetMaxMin ()
 
G4double GetMin ()
 
virtual void ReadFromECHARM (const G4String &, G4double)
 
 ~G4ChannelingECHARM ()
 

Private Attributes

G4double fDistances [3]
 
G4double fMaximum
 
G4double fMinimum
 
G4int fPoints [3]
 
G4PhysicsVectorfVectorEC
 
G4Physics2DVectorfVectorEC2D
 

Detailed Description

Definition at line 34 of file G4ChannelingECHARM.hh.

Constructor & Destructor Documentation

◆ G4ChannelingECHARM()

G4ChannelingECHARM::G4ChannelingECHARM ( const G4String fileName,
G4double  vConversion 
)

Definition at line 32 of file G4ChannelingECHARM.cc.

32 :
33fVectorEC(0),
34fDistances{0.,0.,0.},
35fPoints{0,0,0},
38 fDistances[0] = 0;
39 fDistances[1] = 0;
40 fDistances[2] = 0;
41 fPoints[0] = 0;
42 fPoints[1] = 0;
43 fPoints[2] = 0;
44 fVectorEC2D = 0;
45 ReadFromECHARM(fileName,vConversion);
46}
G4PhysicsVector * fVectorEC
G4Physics2DVector * fVectorEC2D
virtual void ReadFromECHARM(const G4String &, G4double)
#define DBL_MAX
Definition: templates.hh:62

References fDistances, fPoints, fVectorEC2D, and ReadFromECHARM().

◆ ~G4ChannelingECHARM()

G4ChannelingECHARM::~G4ChannelingECHARM ( )

Definition at line 49 of file G4ChannelingECHARM.cc.

49 {
50 delete(fVectorEC);
51 delete(fVectorEC2D);
52}

References fVectorEC, and fVectorEC2D.

Member Function Documentation

◆ GetEC()

G4double G4ChannelingECHARM::GetEC ( G4ThreeVector vPosition)

Definition at line 56 of file G4ChannelingECHARM.cc.

56 {
57 G4double vX = vPosition.x();
58 if (vX < 0.0) {
59 vX += ((int( - vX / fDistances[0]) + 1.0 ) * fDistances[0]);
60 }
61 else if( vX > fDistances[0] ){
62 vX -= ( int( vX / fDistances[0]) * fDistances[0] );
63 }
64 if(fPoints[1]==1){
65 return fVectorEC->Value(vX);
66 }
67 else{
68 G4double vY = vPosition.y();
69 if (vY < 0.0) {
70 vY += ((int( - vY / fDistances[1]) + 1.0 ) * fDistances[1]);
71 }
72 else if( vY > fDistances[1] ){
73 vY -= ( int( vY / fDistances[1]) * fDistances[1] );
74 }
75 return fVectorEC2D->Value((vX),(vY));
76 }
77}
double G4double
Definition: G4Types.hh:83
double x() const
double y() const
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
G4double Value(const G4double energy, std::size_t &lastidx) const

References fDistances, fPoints, fVectorEC, fVectorEC2D, G4PhysicsVector::Value(), G4Physics2DVector::Value(), CLHEP::Hep3Vector::x(), and CLHEP::Hep3Vector::y().

Referenced by G4Channeling::UpdateParameters().

◆ GetIntSp()

G4double G4ChannelingECHARM::GetIntSp ( G4int  index)
inline

Definition at line 49 of file G4ChannelingECHARM.hh.

49{return fDistances[index];};

References fDistances.

Referenced by G4Channeling::UpdateParameters().

◆ GetMax()

G4double G4ChannelingECHARM::GetMax ( )
inline

Definition at line 45 of file G4ChannelingECHARM.hh.

45{return fMaximum;};

References fMaximum.

◆ GetMaxMin()

G4double G4ChannelingECHARM::GetMaxMin ( )
inline

Definition at line 47 of file G4ChannelingECHARM.hh.

47{return std::fabs(fMaximum-fMinimum);};

References fMaximum, and fMinimum.

◆ GetMin()

G4double G4ChannelingECHARM::GetMin ( )
inline

Definition at line 46 of file G4ChannelingECHARM.hh.

46{return fMinimum;};

References fMinimum.

◆ ReadFromECHARM()

void G4ChannelingECHARM::ReadFromECHARM ( const G4String filename,
G4double  vConversion 
)
virtual

Definition at line 81 of file G4ChannelingECHARM.cc.

82 {
83 std::ifstream vFileIn;
84 vFileIn.open(filename);
85
86 vFileIn >> fPoints[0] >> fPoints[1] >> fPoints[2];
87 vFileIn >> fDistances[0] >> fDistances[1] >> fDistances[2];
88
94
95 if(fPoints[1]<1){
97 ed << "No Points not found !" << G4endl;
98 G4Exception("G4ChannelingECHARM::ReadFromECHARM(...)",
99 "G4ChannelingECHARM",
101 ed);
102 return;
103 }
104 else if(fPoints[1]==1){
106 }
107 else{
109 }
110 G4double stepX = fDistances[0]/fPoints[0];
111 G4double stepY = fDistances[1]/fPoints[1];
112 for(G4int i1=0;i1<fPoints[1]; i1++){
113 if(fPoints[1]>1){
114 fVectorEC2D->PutY(i1,i1*stepY);
115 }
116 for(G4int i0=0;i0<fPoints[0]; i0++){
117 double vTempX;
118 vFileIn >> vTempX;
119
120 vTempX *= vConversion;
121 if(vTempX > fMaximum) {fMaximum = vTempX;}
122 if(vTempX < fMinimum) {fMinimum = vTempX;}
123 if(fPoints[1]==1){
124 fVectorEC->PutValue(i0,vTempX);
125 }
126 else{
127 fVectorEC2D->PutValue(i0,i1,vTempX);
128 fVectorEC2D->PutX(i0,i0*stepX);
129 }
130 }
131 }
132 G4cout << "G4ChannelingECHARM::ReadFromECHARM() - " << vConversion << " " << fPoints[0] << " " << fDistances[0] << " " << fPoints[1] << " " << fDistances[1] << " " << fMinimum << " " << fMaximum << G4endl;
133
134 vFileIn.close();
135
136}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void PutY(std::size_t idy, G4double value)
void PutValue(std::size_t idx, std::size_t idy, G4double value)
void PutX(std::size_t idx, G4double value)
void PutValue(const std::size_t index, const G4double value)
static constexpr double meter
Definition: SystemOfUnits.h:71

References DBL_MAX, FatalException, fDistances, fMaximum, fMinimum, fPoints, fVectorEC, fVectorEC2D, G4cout, G4endl, G4Exception(), CLHEP::meter, G4PhysicsVector::PutValue(), G4Physics2DVector::PutValue(), G4Physics2DVector::PutX(), and G4Physics2DVector::PutY().

Referenced by G4ChannelingECHARM().

Field Documentation

◆ fDistances

G4double G4ChannelingECHARM::fDistances[3]
private

Definition at line 55 of file G4ChannelingECHARM.hh.

Referenced by G4ChannelingECHARM(), GetEC(), GetIntSp(), and ReadFromECHARM().

◆ fMaximum

G4double G4ChannelingECHARM::fMaximum
private

Definition at line 58 of file G4ChannelingECHARM.hh.

Referenced by GetMax(), GetMaxMin(), and ReadFromECHARM().

◆ fMinimum

G4double G4ChannelingECHARM::fMinimum
private

Definition at line 59 of file G4ChannelingECHARM.hh.

Referenced by GetMaxMin(), GetMin(), and ReadFromECHARM().

◆ fPoints

G4int G4ChannelingECHARM::fPoints[3]
private

Definition at line 56 of file G4ChannelingECHARM.hh.

Referenced by G4ChannelingECHARM(), GetEC(), and ReadFromECHARM().

◆ fVectorEC

G4PhysicsVector* G4ChannelingECHARM::fVectorEC
private

Definition at line 52 of file G4ChannelingECHARM.hh.

Referenced by GetEC(), ReadFromECHARM(), and ~G4ChannelingECHARM().

◆ fVectorEC2D

G4Physics2DVector* G4ChannelingECHARM::fVectorEC2D
private

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