Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/management/include/G4VITRestProcess.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 //
 27 /// \brief Identical to G4VRestProcess with dependency from G4VITProcess
 28 //
 29 // WARNING : This class is released as a prototype.
 30 // It might strongly evolve or even disapear in the next releases.
 31 //
 32 // Author: Mathieu Karamitros
 33 
 34 // The code is developed in the framework of the ESA AO7146
 35 //
 36 // We would be very happy hearing from you, send us your feedback! :)
 37 //
 38 // In order for Geant4-DNA to be maintained and still open-source,
 39 // article citations are crucial. 
 40 // If you use Geant4-DNA chemistry and you publish papers about your software, 
 41 // in addition to the general paper on Geant4-DNA:
 42 //
 43 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
 44 //
 45 // we would be very happy if you could please also cite the following
 46 // reference papers on chemistry:
 47 //
 48 // J. Comput. Phys. 274 (2014) 841-882
 49 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508 
 50 
 51 #ifndef G4VITRestProcess_h
 52 #define G4VITRestProcess_h 1
 53 
 54 #include <CLHEP/Units/SystemOfUnits.h>
 55 
 56 #include "G4VITProcess.hh"
 57 
 58 /**
 59  * Identical to G4VRestProcess with dependency from G4VITProcess
 60  */
 61 
 62 class G4VITRestProcess : public G4VITProcess
 63 {
 64   //  Abstract class which defines the public behavior of
 65   //  physics interactions at rest.
 66 
 67 public:
 68   G4VITRestProcess(const G4String&, G4ProcessType aType = fNotDefined);
 69   G4VITRestProcess(const G4VITRestProcess&);
 70 
 71   ~G4VITRestProcess() override;
 72 
 73 public:
 74   //  with description
 75   G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
 76                                                       G4ForceCondition* condition) override;
 77 
 78   G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override;
 79 
 80   //  no operation in  PostStepDoIt and  AlongStepDoIt
 81   G4double AlongStepGetPhysicalInteractionLength(const G4Track&,
 82                                                          G4double,
 83                                                          G4double,
 84                                                          G4double&,
 85                                                          G4GPILSelection*) override
 86   {
 87     return -1.0;
 88   }
 89 
 90   G4double PostStepGetPhysicalInteractionLength(const G4Track&,
 91                                                         G4double,
 92                                                         G4ForceCondition*) override
 93   {
 94     return -1.0;
 95   }
 96 
 97   //  no operation in  PostStepDoIt and  AlongStepDoIt
 98   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override
 99   {
100     return nullptr;
101   }
102 
103   G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override
104   {
105     return nullptr;
106   }
107 
108 protected:
109   //  with description
110 
111   virtual G4double GetMeanLifeTime(const G4Track& aTrack,
112                                    G4ForceCondition* condition)=0;
113   //  Calculates the mean life-time (i.e. for decays) of the
114   //  particle at rest due to the occurrence of the given process,
115   //  or converts the probability of interaction (i.e. for
116   //  annihilation) into the life-time of the particle for the
117   //  occurrence of the given process.
118 
119 protected:
120   // hide default constructor and assignment operator as private
121   G4VITRestProcess();
122   G4VITRestProcess & operator=(const G4VITRestProcess &right);
123 };
124 
125 // -----------------------------------------
126 //  inlined function members implementation
127 // -----------------------------------------
128 inline G4double G4VITRestProcess::AtRestGetPhysicalInteractionLength(const G4Track& track,
129                                                                      G4ForceCondition* condition)
130 {
131   // beggining of tracking
132   ResetNumberOfInteractionLengthLeft();
133 
134   // condition is set to "Not Forced"
135   *condition = NotForced;
136 
137   // get mean life time
138   fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
139 
140 #ifdef G4VERBOSE
141   if((fpState->currentInteractionLength < 0.0) || (verboseLevel > 2))
142   {
143     G4cout << "G4VITRestProcess::AtRestGetPhysicalInteractionLength ";
144     G4cout << "[ " << GetProcessName() << "]" << G4endl;
145     track.GetDynamicParticle()->DumpInfo();
146     G4cout << " in Material  " << track.GetMaterial()->GetName() << G4endl;
147     G4cout << "MeanLifeTime = " << fpState->currentInteractionLength / CLHEP::ns
148            << "[ns]" << G4endl;
149   }
150 #endif
151 
152   return (fpState->theNumberOfInteractionLengthLeft)
153       * (fpState->currentInteractionLength);
154 }
155 
156 inline G4VParticleChange* G4VITRestProcess::AtRestDoIt(const G4Track&,
157                                                        const G4Step&)
158 {
159   ClearNumberOfInteractionLengthLeft();
160   ClearInteractionTimeLeft();
161   return pParticleChange;
162 }
163 
164 #endif
165 
166