G4FastTrack.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // $Id:
00030 //---------------------------------------------------------------
00031 //
00032 //  G4FastTrack.hh
00033 //
00034 //  Description:
00035 //    Keeps the current track information and special features
00036 //    for Parameterised Simulation Models.
00037 //
00038 //  History:
00039 //    Oct 97: Verderi && MoraDeFreitas - First Implementation.
00040 //
00041 //---------------------------------------------------------------
00042 
00043 
00044 #ifndef G4FastTrack_h
00045 #define G4FastTrack_h
00046 
00047 #include "G4VSolid.hh"
00048 #include "G4LogicalVolume.hh"
00049 #include "G4Region.hh"
00050 #include "G4AffineTransform.hh"
00051 #include "G4Track.hh"
00052 #include "G4Navigator.hh"
00053 
00054 //---------------------------
00055 // For possible future needs:
00056 //---------------------------
00057 typedef G4Region G4Envelope;
00058 
00059 
00060 //-------------------------------------------
00061 //
00062 //        G4FastTrack class
00063 //
00064 //-------------------------------------------
00065 
00066 // Class Description:
00067 //  The G4FastTrack provides you access to the current G4Track,
00068 // gives simple access to envelope related features (G4Region,
00069 // G4LogicalVolume, G4VSolid, G4AffineTransform references between
00070 // the global and the envelope local coordinates systems) and
00071 // simple access to the position, momentum expressed in the
00072 // envelope coordinate system. Using those quantities and the
00073 // G4VSolid methods, you can for example easily check how far you
00074 // are from the envelope boundary. 
00075 //
00076 
00077 
00078 class G4FastTrack
00079 {
00080 public: // without description
00081   //------------------------
00082   // Constructor/Destructor
00083   //------------------------
00084   // Only one Constructor. By default the envelope can
00085   // be placed n-Times. If the user is sure that it'll be 
00086   // placed just one time, the IsUnique flag should be set 
00087   // TRUE to avoid the G4AffineTransform re-calculations each 
00088   // time we reach the envelope.
00089   G4FastTrack(G4Envelope *anEnvelope,
00090               G4bool IsUnique);
00091   ~G4FastTrack();
00092 
00093   //------------------------------------------------------------
00094   // The fast simulation manager uses the SetCurrentTrack
00095   // method to setup the current G4FastTrack object 
00096   //------------------------------------------------------------
00097   void SetCurrentTrack(const G4Track&, const G4Navigator* a = 0);
00098 
00099   //------------------------------------------------------------
00100   // The fast simulation manager uses the OnTheBoundaryButExiting
00101   // method to test if the particle is leaving the envelope.
00102   //------------------------------------------------------------
00103   G4bool OnTheBoundaryButExiting() const;
00104 
00105   //----------------------------------
00106   // Informations useful to the user :
00107   // General public get functions.
00108   //----------------------------------
00109 
00110 public: // with Description
00111 
00112   const G4Track* GetPrimaryTrack() const;
00113   // Returns the current G4Track.
00114 
00115   G4Envelope* GetEnvelope() const;
00116   // Returns the Envelope G4Region pointer.
00117 
00118   G4LogicalVolume* GetEnvelopeLogicalVolume() const;
00119   // Returns the Envelope G4LogicalVolume pointer.
00120 
00121   G4VPhysicalVolume* GetEnvelopePhysicalVolume() const;
00122   // Returns the Envelope G4VPhysicalVolume pointer.
00123 
00124   G4VSolid* GetEnvelopeSolid() const;
00125   // Returns the Envelope G4VSolid pointer.
00126 
00127   //-----------------------------------
00128   // Primary track informations in the
00129   // Envelope coordinate system.
00130   //-----------------------------------
00131 
00132   G4ThreeVector GetPrimaryTrackLocalPosition() const;
00133   // Returns the particle position in envelope coordinates.
00134 
00135   G4ThreeVector GetPrimaryTrackLocalMomentum() const;
00136   // Returns the particle momentum in envelope coordinates.
00137 
00138   G4ThreeVector GetPrimaryTrackLocalDirection() const;
00139   // Returns the particle direction in envelope coordinates.
00140 
00141   G4ThreeVector GetPrimaryTrackLocalPolarization() const;
00142   // Returns the particle polarization in envelope coordinates.
00143   
00144   //------------------------------------
00145   // 3D transformation of the envelope:
00146   //------------------------------------
00147   // Global -> Local
00148 
00149   const G4AffineTransform* GetAffineTransformation() const;  
00150   // Returns the envelope Global -> Local G4AffineTransform
00151 
00152   // Local -> Global
00153   const G4AffineTransform* GetInverseAffineTransformation() const; 
00154   // Returns the envelope Local -> Global G4AffineTransform
00155 
00156   //-----------------
00157   // Private members
00158   //-----------------
00159 private:
00160 
00161   // Current G4Track pointer
00162   const G4Track* fTrack;
00163 
00164   //------------------------------------------------
00165   // Records the Affine/InverseAffine transformation
00166   // of the envelope.
00167   //------------------------------------------------
00168   void FRecordsAffineTransformation(const G4Navigator*);
00169   G4bool             fAffineTransformationDefined;
00170   G4Envelope*                           fEnvelope;
00171   G4bool                                fIsUnique;
00172   G4LogicalVolume*         fEnvelopeLogicalVolume;
00173   G4VPhysicalVolume*      fEnvelopePhysicalVolume;
00174   G4VSolid*                        fEnvelopeSolid;
00175   G4ThreeVector               fLocalTrackPosition,
00176                               fLocalTrackMomentum, 
00177                              fLocalTrackDirection,
00178                           fLocalTrackPolarization;
00179   G4AffineTransform         fAffineTransformation,
00180                      fInverseAffineTransformation;
00181 };
00182 
00183 
00184 // -----------------
00185 // -- Inline methods
00186 // -----------------
00187 
00188 inline G4Envelope* G4FastTrack::GetEnvelope() const
00189 {
00190   return fEnvelope;
00191 }
00192 
00193 inline G4LogicalVolume* G4FastTrack::GetEnvelopeLogicalVolume() const
00194 {
00195   return fEnvelopeLogicalVolume;
00196 }
00197 
00198 inline G4VPhysicalVolume* G4FastTrack::GetEnvelopePhysicalVolume() const
00199 {
00200   return fEnvelopePhysicalVolume;
00201 }
00202 
00203 inline G4VSolid* G4FastTrack::GetEnvelopeSolid() const
00204 {
00205   return fEnvelopeSolid;
00206 }
00207 
00208 inline const G4Track* G4FastTrack::GetPrimaryTrack() const
00209 {
00210   return fTrack;
00211 }
00212 
00213 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalPosition() const
00214 {
00215   return fLocalTrackPosition;
00216 }
00217 
00218 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalMomentum() const
00219 {
00220   return fLocalTrackMomentum;
00221 }
00222 
00223 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalDirection() const
00224 {
00225   return fLocalTrackDirection;
00226 }
00227 
00228 inline G4ThreeVector G4FastTrack::GetPrimaryTrackLocalPolarization() const
00229 {
00230   return fLocalTrackPolarization;
00231 }
00232 
00233 inline const G4AffineTransform* G4FastTrack::GetAffineTransformation() const
00234 {
00235   return &fAffineTransformation;
00236 }
00237 
00238 inline const G4AffineTransform* G4FastTrack::GetInverseAffineTransformation() const
00239 {
00240   return &fInverseAffineTransformation;
00241 }
00242 
00243 inline G4bool G4FastTrack::OnTheBoundaryButExiting() const 
00244 {
00245   // tests if particle are on the boundary and leaving.
00246   return GetEnvelopeSolid()->
00247     DistanceToOut(GetPrimaryTrackLocalPosition(),
00248                   GetPrimaryTrackLocalDirection())==0.;
00249 }
00250 
00251 #endif

Generated on Mon May 27 17:48:15 2013 for Geant4 by  doxygen 1.4.7