Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/biasing/generic/include/G4ParallelGeometriesLimiterProcess.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // G4ParallelGeometriesLimiterProcess
 27 //
 28 // Description:
 29 //
 30 // Process dedicated to limiting the step on boudaries of
 31 // parallel geometries used for the generic biasing.
 32 //
 33 // Author: Marc Verderi, September 2016.
 34 // --------------------------------------------------------------------
 35 #ifndef G4ParallelGeometriesLimiterProcess_hh
 36 #define G4ParallelGeometriesLimiterProcess_hh 1
 37 
 38 #include "globals.hh"
 39 #include "G4VProcess.hh"
 40 #include "G4Step.hh"
 41 #include "G4Navigator.hh"
 42 #include "G4VPhysicalVolume.hh"
 43 #include "G4ParticleChangeForNothing.hh"
 44 #include "G4FieldTrack.hh"
 45 
 46 class G4PathFinder;
 47 class G4TransportationManager;
 48 
 49 class G4ParallelGeometriesLimiterProcess : public G4VProcess
 50 {
 51   public:
 52 
 53     // --------------------------
 54     // -- Constructor/Destructor:
 55     // --------------------------
 56     G4ParallelGeometriesLimiterProcess(const G4String& processName = "biasLimiter");
 57 
 58     virtual ~G4ParallelGeometriesLimiterProcess() = default;
 59   
 60     // ----------------------------------------------------
 61     // -- Registration / deregistration of parallel worlds.
 62     // -- Access to the list of registered parallel worlds.
 63     // ----------------------------------------------------
 64     void AddParallelWorld(const G4String& parallelWorldName);
 65     void RemoveParallelWorld(const G4String& parallelWorldName);
 66     // -- The list of registered parallel worlds:
 67     const std::vector< G4VPhysicalVolume* >& GetParallelWorlds() const
 68       { return fParallelWorlds; }
 69 
 70     // -- Get the parallel world volume index in the list of world volumes handled
 71     // -- by the process. This index can then be used to access current volume and step
 72     // -- limitation status below.
 73     // -- If the world passed is unknown to the process, -1 is returned.
 74     G4int GetParallelWorldIndex( const G4VPhysicalVolume* parallelWorld ) const;
 75     G4int GetParallelWorldIndex( const G4String& parallelWorldName ) const;
 76 
 77     // ---------------------
 78     // -- Active navigators:
 79     // ---------------------
 80     // -- The list of navigators handled by the process:
 81     const std::vector< G4Navigator* >& GetActiveNavigators() const
 82       { return fParallelWorldNavigators; }
 83     // -- The navigator used for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
 84     // -- Note that no boundary checks are done on the index passed.
 85     const G4Navigator* GetNavigator( G4int worldIndex ) const
 86       { return fParallelWorldNavigators[std::size_t(worldIndex)]; }
 87 
 88     // ---------------------------------------------------
 89     // -- Previous and current steps geometry information:
 90     // ---------------------------------------------------
 91     // -- The "switch" between the previous and current step is done in the PostStepGPIL
 92     // -- The update on the current step is done:
 93     // --    - in the PostStepGPIL for the volumes
 94     // --    - in the AlongStepGPIL for the step limitations
 95     // --
 96     // -- The list of previous step and current step volumes:
 97     const std::vector< const G4VPhysicalVolume* >& GetCurrentVolumes() const
 98       { return  fCurrentVolumes; }
 99     const std::vector< const G4VPhysicalVolume* >& GetPreviousVolumes() const
100       { return fPreviousVolumes; }
101     // -- The current and previous volume for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
102     // -- Note that no boundary checks are done on the index passed.
103     const G4VPhysicalVolume* GetCurrentVolume( G4int worldIndex ) const
104       { return  fCurrentVolumes[std::size_t(worldIndex)]; }
105     const G4VPhysicalVolume* GetPreviousVolume( G4int worldIndex ) const
106       { return fPreviousVolumes[std::size_t(worldIndex)]; }
107     // -- Flags telling about step limitation in previous and current step:
108     const std::vector< G4bool >& GetIsLimiting() const
109       { return  fParallelWorldIsLimiting; }
110     const std::vector< G4bool >& GetWasLimiting() const
111       { return fParallelWorldWasLimiting; }
112     // -- The current and previous step limitation status for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
113     // -- Note that no boundary checks are done on the index passed.
114     G4bool GetIsLimiting( G4int worldIndex ) const
115       { return  fParallelWorldIsLimiting[std::size_t(worldIndex)]; }
116     G4bool GetWasLimiting( G4int worldIndex ) const
117       { return fParallelWorldWasLimiting[std::size_t(worldIndex)]; }
118 
119     // --------------------------------------------------------------
120     //                      From process interface
121     // --------------------------------------------------------------
122     // -- Start/End tracking:
123     void StartTracking(G4Track*);
124     void EndTracking();
125 
126     // --------------------
127     // -- PostStep methods:
128     // --------------------
129     // -- PostStepGPIL is used to collect up to date volumes in the parallel geometries:
130     G4double PostStepGetPhysicalInteractionLength(const G4Track&, G4double, G4ForceCondition*);
131     // -- PostStepDoIt is not used (never called):
132     G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step& ) { return nullptr; }
133 
134     // ---------------------------------------------------------------------------
135     // -- Along step used for limiting the step on parallel geometries boundaries:
136     // ---------------------------------------------------------------------------
137     G4double AlongStepGetPhysicalInteractionLength(const G4Track& track,
138                         G4double previousStepSize, G4double currentMinimumStep, 
139                         G4double& proposedSafety, G4GPILSelection* selection);
140     G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step& step);
141   
142     // ---------------------------
143     // -- AtRest methods not used:
144     // ---------------------------
145     G4double AtRestGetPhysicalInteractionLength(const G4Track&, G4ForceCondition*)
146       { return DBL_MAX; }
147     G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&)
148       { return nullptr; }
149 
150     // --
151     virtual void SetProcessManager(const G4ProcessManager*);
152 
153   private:
154 
155     std::vector< G4VPhysicalVolume* > fParallelWorlds;
156     std::vector< G4Navigator* > fParallelWorldNavigators;
157     std::vector< G4int > fParallelWorldNavigatorIndeces;
158     std::vector< G4double > fParallelWorldSafeties;
159     std::vector< G4bool > fParallelWorldIsLimiting;
160     std::vector< G4bool > fParallelWorldWasLimiting;
161     std::vector< const G4VPhysicalVolume* > fCurrentVolumes;
162     std::vector< const G4VPhysicalVolume* > fPreviousVolumes;
163     G4double fParallelWorldSafety = 0.0;
164     G4bool fIsTrackingTime = false;
165     G4FieldTrack fFieldTrack = '0';
166     G4ParticleChangeForNothing fDummyParticleChange;
167     G4PathFinder* fPathFinder = nullptr;
168     G4TransportationManager* fTransportationManager = nullptr;
169 };
170 
171 #endif
172