Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:16

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 //  G4FastTrack.hh
0031 //
0032 //  Description:
0033 //    Keeps the current track information and special features
0034 //    for Parameterised Simulation Models.
0035 //
0036 //  History:
0037 //    Oct 97: Verderi && MoraDeFreitas - First Implementation.
0038 //
0039 //---------------------------------------------------------------
0040 
0041 #ifndef G4FastTrack_h
0042 #define G4FastTrack_h
0043 
0044 #include "G4AffineTransform.hh"
0045 #include "G4LogicalVolume.hh"
0046 #include "G4Navigator.hh"
0047 #include "G4Region.hh"
0048 #include "G4Track.hh"
0049 #include "G4VSolid.hh"
0050 
0051 //---------------------------
0052 // For possible future needs:
0053 //---------------------------
0054 using G4Envelope = G4Region;
0055 
0056 //-------------------------------------------
0057 //
0058 //        G4FastTrack class
0059 //
0060 //-------------------------------------------
0061 
0062 // Class Description:
0063 //  The G4FastTrack provides you access to the current G4Track,
0064 // gives simple access to envelope related features (G4Region,
0065 // G4LogicalVolume, G4VSolid, G4AffineTransform references between
0066 // the global and the envelope local coordinates systems) and
0067 // simple access to the position, momentum expressed in the
0068 // envelope coordinate system. Using those quantities and the
0069 // G4VSolid methods, you can for example easily check how far you
0070 // are from the envelope boundary.
0071 //
0072 
0073 class G4FastTrack
0074 {
0075   public:  // without description
0076     //------------------------
0077     // Constructor/Destructor
0078     //------------------------
0079     // Only one Constructor. By default the envelope can
0080     // be placed n-Times. If the user is sure that it'll be
0081     // placed just one time, the IsUnique flag should be set
0082     // TRUE to avoid the G4AffineTransform re-calculations each
0083     // time we reach the envelope.
0084     G4FastTrack(G4Envelope* anEnvelope, G4bool IsUnique);
0085     ~G4FastTrack() = default;
0086 
0087     //------------------------------------------------------------
0088     // The fast simulation manager uses the SetCurrentTrack
0089     // method to setup the current G4FastTrack object
0090     //------------------------------------------------------------
0091     void SetCurrentTrack(const G4Track&, const G4Navigator* a = nullptr);
0092 
0093     //------------------------------------------------------------
0094     // The fast simulation manager uses the OnTheBoundaryButExiting
0095     // method to test if the particle is leaving the envelope.
0096     //------------------------------------------------------------
0097     G4bool OnTheBoundaryButExiting() const;
0098 
0099     //----------------------------------
0100     // Informations useful to the user :
0101     // General public get functions.
0102     //----------------------------------
0103 
0104     // Returns the current G4Track.
0105     const G4Track* GetPrimaryTrack() const;
0106 
0107     // Returns the Envelope G4Region pointer.
0108     G4Envelope* GetEnvelope() const;
0109 
0110     // Returns the Envelope G4LogicalVolume pointer.
0111     G4LogicalVolume* GetEnvelopeLogicalVolume() const;
0112 
0113     // Returns the Envelope G4VPhysicalVolume pointer.
0114     G4VPhysicalVolume* GetEnvelopePhysicalVolume() const;
0115 
0116     // Returns the Envelope G4VSolid pointer.
0117     G4VSolid* GetEnvelopeSolid() const;
0118 
0119     //-----------------------------------
0120     // Primary track informations in the
0121     // Envelope coordinate system.
0122     //-----------------------------------
0123 
0124     // Returns the particle position in envelope coordinates.
0125     G4ThreeVector GetPrimaryTrackLocalPosition() const;
0126 
0127     // Returns the particle momentum in envelope coordinates.
0128     G4ThreeVector GetPrimaryTrackLocalMomentum() const;
0129 
0130     // Returns the particle direction in envelope coordinates.
0131     G4ThreeVector GetPrimaryTrackLocalDirection() const;
0132 
0133     // Returns the particle polarization in envelope coordinates.
0134     G4ThreeVector GetPrimaryTrackLocalPolarization() const;
0135 
0136     //------------------------------------
0137     // 3D transformation of the envelope:
0138     //------------------------------------
0139 
0140     // Returns the envelope Global -> Local G4AffineTransform
0141     const G4AffineTransform* GetAffineTransformation() const;
0142 
0143     // Returns the envelope Local -> Global G4AffineTransform
0144     const G4AffineTransform* GetInverseAffineTransformation() const;
0145 
0146   private:
0147     //-----------------
0148     // Private members
0149     //-----------------
0150     // Current G4Track pointer
0151     const G4Track* fTrack{nullptr};
0152 
0153     //------------------------------------------------
0154     // Records the Affine/InverseAffine transformation
0155     // of the envelope.
0156     //------------------------------------------------
0157     void FRecordsAffineTransformation(const G4Navigator*);
0158     G4bool fAffineTransformationDefined{false};
0159     G4Envelope* fEnvelope;
0160     G4bool fIsUnique;
0161     G4LogicalVolume* fEnvelopeLogicalVolume{nullptr};
0162     G4VPhysicalVolume* fEnvelopePhysicalVolume{nullptr};
0163     G4VSolid* fEnvelopeSolid{nullptr};
0164     G4ThreeVector fLocalTrackPosition, fLocalTrackMomentum, fLocalTrackDirection,
0165       fLocalTrackPolarization;
0166     G4AffineTransform fAffineTransformation, fInverseAffineTransformation;
0167 };
0168 
0169 // -----------------
0170 // -- Inline methods
0171 // -----------------
0172 
0173 inline G4Envelope* G4FastTrack::GetEnvelope() const
0174 {
0175   return fEnvelope;
0176 }
0177 
0178 inline G4LogicalVolume* G4FastTrack::GetEnvelopeLogicalVolume() const
0179 {
0180   return fEnvelopeLogicalVolume;
0181 }
0182 
0183 inline G4VPhysicalVolume* G4FastTrack::GetEnvelopePhysicalVolume() const
0184 {
0185   return fEnvelopePhysicalVolume;
0186 }
0187 
0188 inline G4VSolid* G4FastTrack::GetEnvelopeSolid() const
0189 {
0190   return fEnvelopeSolid;
0191 }
0192 
0193 inline const G4Track* G4FastTrack::GetPrimaryTrack() const
0194 {
0195   return fTrack;
0196 }
0197 
0198 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalPosition() const
0199 {
0200   return fLocalTrackPosition;
0201 }
0202 
0203 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalMomentum() const
0204 {
0205   return fLocalTrackMomentum;
0206 }
0207 
0208 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalDirection() const
0209 {
0210   return fLocalTrackDirection;
0211 }
0212 
0213 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalPolarization() const
0214 {
0215   return fLocalTrackPolarization;
0216 }
0217 
0218 inline const G4AffineTransform* G4FastTrack::GetAffineTransformation() const
0219 {
0220   return &fAffineTransformation;
0221 }
0222 
0223 inline const G4AffineTransform* G4FastTrack::GetInverseAffineTransformation() const
0224 {
0225   return &fInverseAffineTransformation;
0226 }
0227 
0228 inline G4bool G4FastTrack::OnTheBoundaryButExiting() const
0229 {
0230   // tests if particle are on the boundary and leaving.
0231   return GetEnvelopeSolid()->DistanceToOut(GetPrimaryTrackLocalPosition(),
0232                                            GetPrimaryTrackLocalDirection())
0233          == 0.;
0234 }
0235 
0236 #endif