Geant4 Cross Reference |
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 // G4BFieldIntegrationDriver implementation 27 // 28 // Specialized integration driver for pure magnetic field 29 // 30 // Author: D.Sorokin 31 // -------------------------------------------------------------------- 32 33 #include "G4BFieldIntegrationDriver.hh" 34 35 #include "G4FieldTrack.hh" 36 #include "G4FieldUtils.hh" 37 #include "G4Exception.hh" 38 #include "G4PhysicalConstants.hh" 39 #include "G4SystemOfUnits.hh" 40 #include "templates.hh" 41 42 43 namespace { 44 45 G4Mag_EqRhs* toMagneticEquation(G4EquationOfMotion* equation) 46 { 47 auto e = dynamic_cast<G4Mag_EqRhs*>(equation); 48 49 if (e == nullptr) 50 { 51 G4Exception("G4BFieldIntegrationDriver::G4BFieldIntegrationDriver", 52 "GeomField0003", FatalErrorInArgument, 53 "Works only with G4Mag_EqRhs"); 54 } 55 56 return e; 57 } 58 59 } // namespace 60 61 62 G4BFieldIntegrationDriver::G4BFieldIntegrationDriver( 63 std::unique_ptr<G4VIntegrationDriver> smallStepDriver, 64 std::unique_ptr<G4VIntegrationDriver> largeStepDriver) 65 : fSmallStepDriver(std::move(smallStepDriver)), 66 fLargeStepDriver(std::move(largeStepDriver)), 67 fCurrDriver(fSmallStepDriver.get()), 68 fEquation(toMagneticEquation(fCurrDriver->GetEquationOfMotion())) 69 { 70 if (fSmallStepDriver->GetEquationOfMotion() 71 != fLargeStepDriver->GetEquationOfMotion()) 72 { 73 G4Exception("G4BFieldIntegrationDriver Constructor:", 74 "GeomField1001", FatalException, "different EoM"); 75 } 76 } 77 78 G4double G4BFieldIntegrationDriver::AdvanceChordLimited(G4FieldTrack& yCurrent, 79 G4double stepMax, 80 G4double epsStep, 81 G4double chordDistance) 82 { 83 const G4double radius = CurvatureRadius(yCurrent); 84 85 G4VIntegrationDriver* driver = nullptr; 86 if (chordDistance < 2 * radius) 87 { 88 stepMax = std::min(stepMax, twopi * radius); 89 driver = fSmallStepDriver.get(); 90 ++fSmallDriverSteps; 91 } else 92 { 93 driver = fLargeStepDriver.get(); 94 ++fLargeDriverSteps; 95 } 96 97 if (driver != fCurrDriver) 98 { 99 driver->OnComputeStep(); 100 } 101 102 fCurrDriver = driver; 103 104 return fCurrDriver->AdvanceChordLimited(yCurrent, stepMax, 105 epsStep, chordDistance); 106 } 107 108 void 109 G4BFieldIntegrationDriver::SetEquationOfMotion(G4EquationOfMotion* equation) 110 { 111 fEquation = toMagneticEquation(equation); 112 fSmallStepDriver->SetEquationOfMotion(equation); 113 fLargeStepDriver->SetEquationOfMotion(equation); 114 } 115 116 G4double 117 G4BFieldIntegrationDriver::CurvatureRadius(const G4FieldTrack& track) const 118 { 119 G4double field[G4Field::MAX_NUMBER_OF_COMPONENTS]; 120 121 GetFieldValue(track, field); 122 123 const G4double Bmag2 = field[0] * field[0] 124 + field[1] * field[1] 125 + field[2] * field[2] ; 126 if (Bmag2 == 0.0 ) 127 { 128 return DBL_MAX; 129 } 130 131 const G4double momentum2 = track.GetMomentum().mag2(); 132 const G4double fCof_inv = eplus / std::abs(fEquation->FCof()); 133 134 return std::sqrt(momentum2 / Bmag2) * fCof_inv; 135 } 136 137 void 138 G4BFieldIntegrationDriver::GetFieldValue(const G4FieldTrack& track, 139 G4double Field[] ) const 140 { 141 G4ThreeVector pos= track.GetPosition(); 142 G4double positionTime[4]= { pos.x(), pos.y(), pos.z(), 143 track.GetLabTimeOfFlight() } ; 144 145 fEquation->GetFieldValue(positionTime, Field); 146 } 147 148 void G4BFieldIntegrationDriver::PrintStatistics() const 149 { 150 const auto totSteps = fSmallDriverSteps + fLargeDriverSteps; 151 const auto toFraction = [&](double value) { return value / totSteps * 100; }; 152 153 G4cout << "============= G4BFieldIntegrationDriver statistics ===========\n" 154 << "total steps " << totSteps << " " 155 << "smallDriverSteps " << toFraction(fSmallDriverSteps) << " " 156 << "largeDriverSteps " << toFraction(fLargeDriverSteps) << "\n" 157 << "======================================\n"; 158 } 159