Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-15 08:59:08

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 // G4ParallelGeometriesLimiterProcess
0027 //
0028 // Description:
0029 //
0030 // Process dedicated to limiting the step on boudaries of
0031 // parallel geometries used for the generic biasing.
0032 //
0033 // Author: Marc Verderi, September 2016.
0034 // --------------------------------------------------------------------
0035 #ifndef G4ParallelGeometriesLimiterProcess_hh
0036 #define G4ParallelGeometriesLimiterProcess_hh 1
0037 
0038 #include "globals.hh"
0039 #include "G4VProcess.hh"
0040 #include "G4Step.hh"
0041 #include "G4Navigator.hh"
0042 #include "G4VPhysicalVolume.hh"
0043 #include "G4ParticleChangeForNothing.hh"
0044 #include "G4FieldTrack.hh"
0045 
0046 class G4PathFinder;
0047 class G4TransportationManager;
0048 
0049 class G4ParallelGeometriesLimiterProcess : public G4VProcess
0050 {
0051   public:
0052 
0053     // --------------------------
0054     // -- Constructor/Destructor:
0055     // --------------------------
0056     G4ParallelGeometriesLimiterProcess(const G4String& processName = "biasLimiter");
0057 
0058     virtual ~G4ParallelGeometriesLimiterProcess() = default;
0059   
0060     // ----------------------------------------------------
0061     // -- Registration / deregistration of parallel worlds.
0062     // -- Access to the list of registered parallel worlds.
0063     // ----------------------------------------------------
0064     void AddParallelWorld(const G4String& parallelWorldName);
0065     void RemoveParallelWorld(const G4String& parallelWorldName);
0066     // -- The list of registered parallel worlds:
0067     const std::vector< G4VPhysicalVolume* >& GetParallelWorlds() const
0068       { return fParallelWorlds; }
0069 
0070     // -- Get the parallel world volume index in the list of world volumes handled
0071     // -- by the process. This index can then be used to access current volume and step
0072     // -- limitation status below.
0073     // -- If the world passed is unknown to the process, -1 is returned.
0074     G4int GetParallelWorldIndex( const G4VPhysicalVolume* parallelWorld ) const;
0075     G4int GetParallelWorldIndex( const G4String& parallelWorldName ) const;
0076 
0077     // ---------------------
0078     // -- Active navigators:
0079     // ---------------------
0080     // -- The list of navigators handled by the process:
0081     const std::vector< G4Navigator* >& GetActiveNavigators() const
0082       { return fParallelWorldNavigators; }
0083     // -- The navigator used for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
0084     // -- Note that no boundary checks are done on the index passed.
0085     const G4Navigator* GetNavigator( G4int worldIndex ) const
0086       { return fParallelWorldNavigators[std::size_t(worldIndex)]; }
0087 
0088     // ---------------------------------------------------
0089     // -- Previous and current steps geometry information:
0090     // ---------------------------------------------------
0091     // -- The "switch" between the previous and current step is done in the PostStepGPIL
0092     // -- The update on the current step is done:
0093     // --    - in the PostStepGPIL for the volumes
0094     // --    - in the AlongStepGPIL for the step limitations
0095     // --
0096     // -- The list of previous step and current step volumes:
0097     const std::vector< const G4VPhysicalVolume* >& GetCurrentVolumes() const
0098       { return  fCurrentVolumes; }
0099     const std::vector< const G4VPhysicalVolume* >& GetPreviousVolumes() const
0100       { return fPreviousVolumes; }
0101     // -- The current and previous volume for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
0102     // -- Note that no boundary checks are done on the index passed.
0103     const G4VPhysicalVolume* GetCurrentVolume( G4int worldIndex ) const
0104       { return  fCurrentVolumes[std::size_t(worldIndex)]; }
0105     const G4VPhysicalVolume* GetPreviousVolume( G4int worldIndex ) const
0106       { return fPreviousVolumes[std::size_t(worldIndex)]; }
0107     // -- Flags telling about step limitation in previous and current step:
0108     const std::vector< G4bool >& GetIsLimiting() const
0109       { return  fParallelWorldIsLimiting; }
0110     const std::vector< G4bool >& GetWasLimiting() const
0111       { return fParallelWorldWasLimiting; }
0112     // -- The current and previous step limitation status for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
0113     // -- Note that no boundary checks are done on the index passed.
0114     G4bool GetIsLimiting( G4int worldIndex ) const
0115       { return  fParallelWorldIsLimiting[std::size_t(worldIndex)]; }
0116     G4bool GetWasLimiting( G4int worldIndex ) const
0117       { return fParallelWorldWasLimiting[std::size_t(worldIndex)]; }
0118 
0119     // --------------------------------------------------------------
0120     //                      From process interface
0121     // --------------------------------------------------------------
0122     // -- Start/End tracking:
0123     void StartTracking(G4Track*);
0124     void EndTracking();
0125 
0126     // --------------------
0127     // -- PostStep methods:
0128     // --------------------
0129     // -- PostStepGPIL is used to collect up to date volumes in the parallel geometries:
0130     G4double PostStepGetPhysicalInteractionLength(const G4Track&, G4double, G4ForceCondition*);
0131     // -- PostStepDoIt is not used (never called):
0132     G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step& ) { return nullptr; }
0133 
0134     // ---------------------------------------------------------------------------
0135     // -- Along step used for limiting the step on parallel geometries boundaries:
0136     // ---------------------------------------------------------------------------
0137     G4double AlongStepGetPhysicalInteractionLength(const G4Track& track,
0138                         G4double previousStepSize, G4double currentMinimumStep, 
0139                         G4double& proposedSafety, G4GPILSelection* selection);
0140     G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step& step);
0141   
0142     // ---------------------------
0143     // -- AtRest methods not used:
0144     // ---------------------------
0145     G4double AtRestGetPhysicalInteractionLength(const G4Track&, G4ForceCondition*)
0146       { return DBL_MAX; }
0147     G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&)
0148       { return nullptr; }
0149 
0150     // --
0151     virtual void SetProcessManager(const G4ProcessManager*);
0152 
0153   private:
0154 
0155     std::vector< G4VPhysicalVolume* > fParallelWorlds;
0156     std::vector< G4Navigator* > fParallelWorldNavigators;
0157     std::vector< G4int > fParallelWorldNavigatorIndeces;
0158     std::vector< G4double > fParallelWorldSafeties;
0159     std::vector< G4bool > fParallelWorldIsLimiting;
0160     std::vector< G4bool > fParallelWorldWasLimiting;
0161     std::vector< const G4VPhysicalVolume* > fCurrentVolumes;
0162     std::vector< const G4VPhysicalVolume* > fPreviousVolumes;
0163     G4double fParallelWorldSafety = 0.0;
0164     G4bool fIsTrackingTime = false;
0165     G4FieldTrack fFieldTrack = '0';
0166     G4ParticleChangeForNothing fDummyParticleChange;
0167     G4PathFinder* fPathFinder = nullptr;
0168     G4TransportationManager* fTransportationManager = nullptr;
0169 };
0170 
0171 #endif