Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:24

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 /// \brief Identical to G4VRestProcess with dependency from G4VITProcess
0028 //
0029 // WARNING : This class is released as a prototype.
0030 // It might strongly evolve or even disapear in the next releases.
0031 //
0032 // Author: Mathieu Karamitros
0033 
0034 // The code is developed in the framework of the ESA AO7146
0035 //
0036 // We would be very happy hearing from you, send us your feedback! :)
0037 //
0038 // In order for Geant4-DNA to be maintained and still open-source,
0039 // article citations are crucial. 
0040 // If you use Geant4-DNA chemistry and you publish papers about your software, 
0041 // in addition to the general paper on Geant4-DNA:
0042 //
0043 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
0044 //
0045 // we would be very happy if you could please also cite the following
0046 // reference papers on chemistry:
0047 //
0048 // J. Comput. Phys. 274 (2014) 841-882
0049 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508 
0050 
0051 #ifndef G4VITRestProcess_h
0052 #define G4VITRestProcess_h 1
0053 
0054 #include <CLHEP/Units/SystemOfUnits.h>
0055 
0056 #include "G4VITProcess.hh"
0057 
0058 /**
0059  * Identical to G4VRestProcess with dependency from G4VITProcess
0060  */
0061 
0062 class G4VITRestProcess : public G4VITProcess
0063 {
0064   //  Abstract class which defines the public behavior of
0065   //  physics interactions at rest.
0066 
0067 public:
0068   G4VITRestProcess(const G4String&, G4ProcessType aType = fNotDefined);
0069   G4VITRestProcess(const G4VITRestProcess&);
0070 
0071   ~G4VITRestProcess() override;
0072 
0073 public:
0074   //  with description
0075   G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
0076                                                       G4ForceCondition* condition) override;
0077 
0078   G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override;
0079 
0080   //  no operation in  PostStepDoIt and  AlongStepDoIt
0081   G4double AlongStepGetPhysicalInteractionLength(const G4Track&,
0082                                                          G4double,
0083                                                          G4double,
0084                                                          G4double&,
0085                                                          G4GPILSelection*) override
0086   {
0087     return -1.0;
0088   }
0089 
0090   G4double PostStepGetPhysicalInteractionLength(const G4Track&,
0091                                                         G4double,
0092                                                         G4ForceCondition*) override
0093   {
0094     return -1.0;
0095   }
0096 
0097   //  no operation in  PostStepDoIt and  AlongStepDoIt
0098   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override
0099   {
0100     return nullptr;
0101   }
0102 
0103   G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override
0104   {
0105     return nullptr;
0106   }
0107 
0108 protected:
0109   //  with description
0110 
0111   virtual G4double GetMeanLifeTime(const G4Track& aTrack,
0112                                    G4ForceCondition* condition)=0;
0113   //  Calculates the mean life-time (i.e. for decays) of the
0114   //  particle at rest due to the occurrence of the given process,
0115   //  or converts the probability of interaction (i.e. for
0116   //  annihilation) into the life-time of the particle for the
0117   //  occurrence of the given process.
0118 
0119 protected:
0120   // hide default constructor and assignment operator as private
0121   G4VITRestProcess();
0122   G4VITRestProcess & operator=(const G4VITRestProcess &right);
0123 };
0124 
0125 // -----------------------------------------
0126 //  inlined function members implementation
0127 // -----------------------------------------
0128 inline G4double G4VITRestProcess::AtRestGetPhysicalInteractionLength(const G4Track& track,
0129                                                                      G4ForceCondition* condition)
0130 {
0131   // beggining of tracking
0132   ResetNumberOfInteractionLengthLeft();
0133 
0134   // condition is set to "Not Forced"
0135   *condition = NotForced;
0136 
0137   // get mean life time
0138   fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
0139 
0140 #ifdef G4VERBOSE
0141   if((fpState->currentInteractionLength < 0.0) || (verboseLevel > 2))
0142   {
0143     G4cout << "G4VITRestProcess::AtRestGetPhysicalInteractionLength ";
0144     G4cout << "[ " << GetProcessName() << "]" << G4endl;
0145     track.GetDynamicParticle()->DumpInfo();
0146     G4cout << " in Material  " << track.GetMaterial()->GetName() << G4endl;
0147     G4cout << "MeanLifeTime = " << fpState->currentInteractionLength / CLHEP::ns
0148            << "[ns]" << G4endl;
0149   }
0150 #endif
0151 
0152   return (fpState->theNumberOfInteractionLengthLeft)
0153       * (fpState->currentInteractionLength);
0154 }
0155 
0156 inline G4VParticleChange* G4VITRestProcess::AtRestDoIt(const G4Track&,
0157                                                        const G4Step&)
0158 {
0159   ClearNumberOfInteractionLengthLeft();
0160   ClearInteractionTimeLeft();
0161   return pParticleChange;
0162 }
0163 
0164 #endif
0165