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

G4DriverReporter is an auxiliary utility class to print information from integration drivers. It can be used by different types of drivers. More...

#include <G4DriverReporter.hh>

Static Public Member Functions

static void PrintStatus (const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, unsigned int subStepNo, unsigned int noIntegrationVariables)
static void PrintStatus (const G4FieldTrack &StartFT, const G4FieldTrack &CurrentFT, G4double requestStep, unsigned int subStepNo)
static void PrintStat_Aux (const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)

Detailed Description

G4DriverReporter is an auxiliary utility class to print information from integration drivers. It can be used by different types of drivers.

Definition at line 45 of file G4DriverReporter.hh.

Member Function Documentation

◆ PrintStat_Aux()

void G4DriverReporter::PrintStat_Aux ( const G4FieldTrack & aFieldTrack,
G4double requestStep,
G4double actualStep,
G4int subStepNo,
G4double subStepSize,
G4double dotVelocities )
static

Definition at line 132 of file G4DriverReporter.cc.

138{
139 const G4int noPrecision = 8;
140 const G4int prec7 = noPrecision+2;
141 const G4int prec8 = noPrecision+3;
142 const G4int prec9 = noPrecision+4;
143
144 const G4ThreeVector Position = aFieldTrack.GetPosition();
145 const G4ThreeVector UnitVelocity = aFieldTrack.GetMomentumDir();
146
147 G4long oldprec= G4cout.precision(noPrecision);
148
149 if( subStepNo >= 0)
150 {
151 G4cout << std::setw( 5) << subStepNo << " ";
152 }
153 else
154 {
155 G4cout << std::setw( 5) << "Start" << " ";
156 }
157 G4double curveLen= aFieldTrack.GetCurveLength();
158 G4cout << std::setw( 7) << curveLen;
159 // G4cout.precision(noPrecision);
160 G4cout << std::setw( prec9) << Position.x() << " "
161 << std::setw( prec9) << Position.y() << " "
162 << std::setw( prec9) << Position.z() << " "
163 << std::setw( prec8) << UnitVelocity.x() << " "
164 << std::setw( prec8) << UnitVelocity.y() << " "
165 << std::setw( prec8) << UnitVelocity.z() << " ";
166 G4cout.precision(3);
167 G4double unitMagDif = UnitVelocity.mag2()-1.0;
168 if( std::fabs( unitMagDif ) < 1.0e-15 ) { unitMagDif = 0.0; }
169 G4cout << std::setw( 8) << unitMagDif << " ";
170 G4cout.precision(6);
171 G4cout << std::setw(10) << dotVeloc_StartCurr << " ";
172 G4cout.precision(oldprec);
173 G4cout << std::setw( prec7) << aFieldTrack.GetKineticEnergy();
174 G4cout << std::setw(12) << step_len << " ";
175
176 static G4ThreadLocal G4double oldCurveLength = 0.0;
177 static G4ThreadLocal G4double oldSubStepLength = 0.0;
178 static G4ThreadLocal G4int oldSubStepNo = -1;
179
180 G4double subStep_len = 0.0;
181 if( curveLen > oldCurveLength )
182 {
183 subStep_len= curveLen - oldCurveLength;
184 }
185 else if (subStepNo == oldSubStepNo)
186 {
187 subStep_len= oldSubStepLength;
188 }
189 oldCurveLength= curveLen;
190 oldSubStepLength= subStep_len;
191
192 G4cout << std::setw(12) << subStep_len << " ";
193 G4cout << std::setw(12) << subStepSize << " ";
194 if( requestStep != -1.0 )
195 {
196 G4cout << std::setw( prec9) << requestStep << " ";
197 }
198 else
199 {
200 G4cout << std::setw( prec9) << " InitialStep " << " ";
201 }
202 G4cout << G4endl;
203}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double mag2() const
double y() const
const G4ThreeVector & GetMomentumDir() const
G4double GetCurveLength() const
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
#define G4ThreadLocal
Definition tls.hh:77

Referenced by PrintStatus().

◆ PrintStatus() [1/2]

void G4DriverReporter::PrintStatus ( const G4double * StartArr,
G4double xstart,
const G4double * CurrentArr,
G4double xcurrent,
G4double requestStep,
unsigned int subStepNo,
unsigned int noIntegrationVariables )
static

Definition at line 35 of file G4DriverReporter.cc.

46{
47 G4FieldTrack StartFT(G4ThreeVector(0,0,0),
48 G4ThreeVector(0,0,0), 0., 0., 0., 0. );
49 G4FieldTrack CurrentFT (StartFT);
50
51 StartFT.LoadFromArray( StartArr, noIntegrationVariables);
52 StartFT.SetCurveLength( xstart);
53 CurrentFT.LoadFromArray( CurrentArr, noIntegrationVariables);
54 CurrentFT.SetCurveLength( xcurrent );
55
56 PrintStatus(StartFT, CurrentFT, requestStep, subStepNo );
57}
static void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, unsigned int subStepNo, unsigned int noIntegrationVariables)

Referenced by G4MagInt_Driver::AccurateAdvance(), G4OldMagIntDriver::AccurateAdvance(), and PrintStatus().

◆ PrintStatus() [2/2]

void G4DriverReporter::PrintStatus ( const G4FieldTrack & StartFT,
const G4FieldTrack & CurrentFT,
G4double requestStep,
unsigned int subStepNo )
static

Definition at line 61 of file G4DriverReporter.cc.

65{
66 const G4int noPrecision = 8;
67 const G4int prec7 = noPrecision+2;
68 const G4int prec8 = noPrecision+3;
69 const G4int prec9 = noPrecision+4;
70
71 const G4ThreeVector StartPosition= StartFT.GetPosition();
72 const G4ThreeVector StartUnitVelocity= StartFT.GetMomentumDir();
73 const G4ThreeVector CurrentPosition= CurrentFT.GetPosition();
74 const G4ThreeVector CurrentUnitVelocity= CurrentFT.GetMomentumDir();
75
76 G4int verboseLevel= 2; // fVerboseLevel;
77 G4long oldPrec= G4cout.precision(noPrecision);
78
79 G4double DotStartCurrentVeloc= StartUnitVelocity.dot(CurrentUnitVelocity);
80
81 G4double step_len= CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
82 G4double subStepSize = step_len;
83
84 if( (subStepNo <= 1) || (verboseLevel > 3) )
85 {
86 subStepNo = - subStepNo; // To allow printing banner
87
88 G4cout << "------------------------------------------------------------------"
89 << G4endl;
90 G4cout << std::setw( 6) << " " << std::setw( 25)
91 << " G4DriverReporter: Current Position and Direction" << " "
92 << G4endl;
93 G4cout << std::setw( 5) << "Step#" << " "
94 << std::setw( prec7) << "s-curve" << " "
95 << std::setw( prec9) << "X(mm)" << " "
96 << std::setw( prec9) << "Y(mm)" << " "
97 << std::setw( prec9) << "Z(mm)" << " "
98 << std::setw( prec8) << " N_x " << " "
99 << std::setw( prec8) << " N_y " << " "
100 << std::setw( prec8) << " N_z " << " "
101 << std::setw( 6) << " N^2-1 " << " "
102 << std::setw(10) << " N(0).N " << " "
103 << std::setw( 7) << "KinEner " << " "
104 << std::setw(12) << "Track-l" << " " // Add the Sub-step ??
105 << std::setw(12) << "Step-len" << " "
106 << std::setw(12) << "Step-len" << " "
107 << std::setw( 9) << "ReqStep" << " "
108 << G4endl;
109 }
110
111 G4cout.precision(noPrecision);
112
113 if( (subStepNo <= 0) )
114 {
115 PrintStat_Aux( StartFT, requestStep, 0.,
116 0, 0.0, 1.0);
117 }
118
119 // if( verboseLevel <= 3 )
120 {
121 G4cout.precision(noPrecision);
122 PrintStat_Aux( CurrentFT, requestStep, step_len,
123 subStepNo, subStepSize, DotStartCurrentVeloc );
124 }
125 G4cout << "------------------------------------------------------------------"
126 << G4endl;
127 G4cout.precision(oldPrec);
128}
double dot(const Hep3Vector &) const
static void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)

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