|
||||
File indexing completed on 2025-01-18 09:58:50
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 0028 // 0029 //--------------------------------------------------------------- 0030 // 0031 // G4ParallelGeometriesLimiterProcess.hh 0032 // 0033 // Description: 0034 // Process dedicated to limiting the step on boudaries of 0035 // parallel geometries used for the generic biasing. 0036 // 0037 // History: 0038 // Sep 16: Creation. 0039 // 0040 //--------------------------------------------------------------- 0041 0042 0043 #ifndef G4ParallelGeometriesLimiterProcess_hh 0044 #define G4ParallelGeometriesLimiterProcess_hh 0045 0046 #include "globals.hh" 0047 #include "G4VProcess.hh" 0048 #include "G4Step.hh" 0049 #include "G4Navigator.hh" 0050 #include "G4VPhysicalVolume.hh" 0051 #include "G4ParticleChangeForNothing.hh" 0052 #include "G4FieldTrack.hh" 0053 class G4PathFinder; 0054 class G4TransportationManager; 0055 0056 0057 // Class Description: 0058 // -- G4VProcess limiting the step on parallel geometries used by generic biasing. 0059 0060 0061 class G4ParallelGeometriesLimiterProcess : public G4VProcess 0062 { 0063 public: 0064 // -------------------------- 0065 // -- Constructor/Destructor: 0066 // -------------------------- 0067 G4ParallelGeometriesLimiterProcess(const G4String& processName = "biasLimiter"); 0068 0069 public: 0070 virtual ~G4ParallelGeometriesLimiterProcess() 0071 {} 0072 0073 // ---------------------------------------------------- 0074 // -- Registration / deregistration of parallel worlds. 0075 // -- Access to the list of registered parallel worlds. 0076 // ---------------------------------------------------- 0077 void AddParallelWorld(const G4String& parallelWorldName); 0078 void RemoveParallelWorld(const G4String& parallelWorldName); 0079 // -- The list of registered parallel worlds: 0080 const std::vector< G4VPhysicalVolume* >& GetParallelWorlds() const { return fParallelWorlds; } 0081 0082 // -- Get the parallel world volume index in the list of world volumes handled 0083 // -- by the process. This index can then be used to access current volume and step 0084 // -- limitation status below. 0085 // -- If the world passed is unknown to the process, -1 is returned. 0086 G4int GetParallelWorldIndex( const G4VPhysicalVolume* parallelWorld ) const; 0087 G4int GetParallelWorldIndex( G4String parallelWorldName ) const; 0088 0089 0090 // --------------------- 0091 // -- Active navigators: 0092 // --------------------- 0093 // -- The list of navigators handled by the process: 0094 const std::vector< G4Navigator* >& GetActiveNavigators() const { return fParallelWorldNavigators; } 0095 // -- The navigator used for the passed parallel world index (obtained with GetParallelWorldIndex(...) above) 0096 // -- Note that no boundary checks are done on the index passed. 0097 const G4Navigator* GetNavigator( G4int worldIndex ) const { return fParallelWorldNavigators[size_t(worldIndex)]; } 0098 // --------------------------------------------------- 0099 // -- Previous and current steps geometry information: 0100 // --------------------------------------------------- 0101 // -- The "switch" between the previous and current step is done in the PostStepGPIL 0102 // -- The update on the current step is done: 0103 // -- - in the PostStepGPIL for the volumes 0104 // -- - in the AlongStepGPIL for the step limitations 0105 // -- 0106 // -- The list of previous step and current step volumes: 0107 const std::vector< const G4VPhysicalVolume* >& GetCurrentVolumes() const { return fCurrentVolumes; } 0108 const std::vector< const G4VPhysicalVolume* >& GetPreviousVolumes() const { return fPreviousVolumes; } 0109 // -- The current and previous volume for the passed parallel world index (obtained with GetParallelWorldIndex(...) above) 0110 // -- Note that no boundary checks are done on the index passed. 0111 const G4VPhysicalVolume* GetCurrentVolume( G4int worldIndex ) const { return fCurrentVolumes[size_t(worldIndex)]; } 0112 const G4VPhysicalVolume* GetPreviousVolume( G4int worldIndex ) const { return fPreviousVolumes[size_t(worldIndex)]; } 0113 // -- Flags telling about step limitation in previous and current step: 0114 const std::vector< G4bool >& GetIsLimiting() const { return fParallelWorldIsLimiting; } 0115 const std::vector< G4bool >& GetWasLimiting() const { return fParallelWorldWasLimiting; } 0116 // -- The current and previous step limitation status for the passed parallel world index (obtained with GetParallelWorldIndex(...) above) 0117 // -- Note that no boundary checks are done on the index passed. 0118 G4bool GetIsLimiting( G4int worldIndex ) const { return fParallelWorldIsLimiting[size_t(worldIndex)]; } 0119 G4bool GetWasLimiting( G4int worldIndex ) const { return fParallelWorldWasLimiting[size_t(worldIndex)]; } 0120 0121 0122 // -------------------------------------------------------------- 0123 // From process interface 0124 // -------------------------------------------------------------- 0125 // -- Start/End tracking: 0126 void StartTracking(G4Track*); 0127 void EndTracking(); 0128 0129 // -------------------- 0130 // -- PostStep methods: 0131 // -------------------- 0132 // -- PostStepGPIL is used to collect up to date volumes in the parallel geometries: 0133 G4double PostStepGetPhysicalInteractionLength(const G4Track&, G4double, G4ForceCondition*); 0134 // -- PostStepDoIt is not used (never called): 0135 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step& ) 0136 { return nullptr; } 0137 0138 // --------------------------------------------------------------------------- 0139 // -- Along step used for limiting the step on parallel geometries boundaries: 0140 // --------------------------------------------------------------------------- 0141 G4double AlongStepGetPhysicalInteractionLength(const G4Track& track, 0142 G4double previousStepSize, 0143 G4double currentMinimumStep, 0144 G4double& proposedSafety, 0145 G4GPILSelection* selection); 0146 G4VParticleChange* AlongStepDoIt(const G4Track& track, 0147 const G4Step& step); 0148 0149 // --------------------------- 0150 // -- AtRest methods not used: 0151 // --------------------------- 0152 G4double AtRestGetPhysicalInteractionLength(const G4Track&, 0153 G4ForceCondition*) 0154 { return DBL_MAX; } 0155 0156 G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) 0157 { return nullptr; } 0158 0159 0160 // -- 0161 virtual void SetProcessManager(const G4ProcessManager*); 0162 0163 0164 0165 private: 0166 std::vector< G4VPhysicalVolume* > fParallelWorlds; 0167 std::vector< G4Navigator* > fParallelWorldNavigators; 0168 std::vector< G4int > fParallelWorldNavigatorIndeces; 0169 std::vector< G4double > fParallelWorldSafeties; 0170 std::vector< G4bool > fParallelWorldIsLimiting; 0171 std::vector< G4bool > fParallelWorldWasLimiting; 0172 std::vector< const G4VPhysicalVolume* > fCurrentVolumes; 0173 std::vector< const G4VPhysicalVolume* > fPreviousVolumes; 0174 G4double fParallelWorldSafety; 0175 G4bool fIsTrackingTime; 0176 G4FieldTrack fFieldTrack; 0177 G4ParticleChangeForNothing fDummyParticleChange; 0178 G4PathFinder* fPathFinder; 0179 G4TransportationManager* fTransportationManager; 0180 }; 0181 0182 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |