Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/polarisation/include/G4PolarizedAnnihilationModel.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 //
 28 // Geant4 Class header file
 29 //
 30 // File name:     G4PolarizedAnnihilationModel
 31 //
 32 // Author:        Andreas Schaelicke and Pavel Starovoitov
 33 //
 34 // Class Description:
 35 //   Implementation of polarized gamma Annihilation scattering on free electron
 36 //
 37 // -------------------------------------------------------------------
 38 
 39 #ifndef G4PolarizedAnnihilationModel_h
 40 #define G4PolarizedAnnihilationModel_h 1
 41 
 42 #include "globals.hh"
 43 #include "G4eeToTwoGammaModel.hh"
 44 #include "G4ThreeVector.hh"
 45 #include "G4StokesVector.hh"
 46 
 47 class G4DynamicParticle;
 48 class G4MaterialCutsCouple;
 49 class G4ParticleChangeForGamma;
 50 class G4ParticleDefinition;
 51 class G4PolarizedAnnihilationXS;
 52 
 53 class G4PolarizedAnnihilationModel : public G4eeToTwoGammaModel
 54 {
 55  public:
 56   explicit G4PolarizedAnnihilationModel(
 57     const G4ParticleDefinition* p = nullptr,
 58     const G4String& nam           = "Polarized-Annihilation");
 59 
 60   virtual ~G4PolarizedAnnihilationModel() override;
 61 
 62   virtual void Initialise(const G4ParticleDefinition*,
 63                           const G4DataVector&) override;
 64 
 65   virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy) override;
 66 
 67   void ComputeAsymmetriesPerElectron(G4double gammaEnergy, G4double& valueX,
 68                                      G4double& valueA, G4double& valueT);
 69 
 70   virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 71                                  const G4MaterialCutsCouple*,
 72                                  const G4DynamicParticle*, G4double tmin,
 73                                  G4double maxEnergy) final;
 74 
 75   // polarized routines
 76   inline void SetTargetPolarization(const G4ThreeVector& pTarget);
 77   inline void SetBeamPolarization(const G4ThreeVector& pBeam);
 78   inline const G4ThreeVector& GetTargetPolarization() const;
 79   inline const G4ThreeVector& GetBeamPolarization() const;
 80   inline const G4ThreeVector& GetFinalGamma1Polarization() const;
 81   inline const G4ThreeVector& GetFinalGamma2Polarization() const;
 82 
 83   G4PolarizedAnnihilationModel& operator                            =(
 84     const G4PolarizedAnnihilationModel& right) = delete;
 85   G4PolarizedAnnihilationModel(const G4PolarizedAnnihilationModel&) = delete;
 86 
 87  private:
 88   G4PolarizedAnnihilationXS* fCrossSectionCalculator;
 89   G4ParticleChangeForGamma* fParticleChange;
 90 
 91   // incoming
 92   G4StokesVector fBeamPolarization;    // positron
 93   G4StokesVector fTargetPolarization;  // electron
 94   // outgoing
 95   G4StokesVector fFinalGamma1Polarization;
 96   G4StokesVector fFinalGamma2Polarization;
 97 
 98   G4int fVerboseLevel;
 99 };
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102 
103 inline void G4PolarizedAnnihilationModel::SetTargetPolarization(
104   const G4ThreeVector& pTarget)
105 {
106   fTargetPolarization = G4StokesVector(pTarget);
107 }
108 inline void G4PolarizedAnnihilationModel::SetBeamPolarization(
109   const G4ThreeVector& pBeam)
110 {
111   fBeamPolarization = G4StokesVector(pBeam);
112 }
113 inline const G4ThreeVector&
114 G4PolarizedAnnihilationModel::GetTargetPolarization() const
115 {
116   return fTargetPolarization;
117 }
118 inline const G4ThreeVector& G4PolarizedAnnihilationModel::GetBeamPolarization()
119   const
120 {
121   return fBeamPolarization;
122 }
123 inline const G4ThreeVector&
124 G4PolarizedAnnihilationModel::GetFinalGamma1Polarization() const
125 {
126   return fFinalGamma1Polarization;
127 }
128 inline const G4ThreeVector&
129 G4PolarizedAnnihilationModel::GetFinalGamma2Polarization() const
130 {
131   return fFinalGamma2Polarization;
132 }
133 
134 #endif
135