Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/Geant4/G4PolarizedAnnihilationXS.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

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 // Geant4 Class file
0029 //
0030 // File name:     G4PolarizedAnnihilationXS
0031 //
0032 // Author:        Andreas Schaelicke and Pavel Starovoitov
0033 //
0034 // Class Description:
0035 //   * calculates the differential cross section (ME squared,
0036 //     without phase space) incoming positron (along positive z direction)
0037 //     annihilations with an electron at rest
0038 //   * phi denotes the angle between the scattering plane and
0039 //     X axis of incoming partice reference frame (PRF)
0040 
0041 #ifndef G4PolarizedAnnihilationXS_h
0042 #define G4PolarizedAnnihilationXS_h 1
0043 
0044 #include "G4StokesVector.hh"
0045 #include "G4VPolarizedXS.hh"
0046 
0047 class G4PolarizedAnnihilationXS : public G4VPolarizedXS
0048 {
0049  public:
0050   G4PolarizedAnnihilationXS();
0051   virtual ~G4PolarizedAnnihilationXS() override;
0052 
0053   virtual void Initialize(G4double eps, G4double gamma, G4double phi,
0054                           const G4StokesVector& p0, const G4StokesVector& p1,
0055                           G4int flag = 0) override;
0056 
0057   G4double DiceEpsilon();
0058 
0059   virtual G4double XSection(const G4StokesVector& pol2,
0060                             const G4StokesVector& pol3) override;
0061 
0062   virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y,
0063                                  const G4StokesVector& pol0,
0064                                  const G4StokesVector& pol1) override;
0065 
0066   // return expected mean polarisation
0067   virtual G4StokesVector GetPol2() override;
0068   virtual G4StokesVector GetPol3() override;
0069 
0070   // minimal energy fraction in TotalXSection
0071   virtual G4double GetXmin(G4double y) override;
0072   // maximal energy fraction in TotalXSection
0073   virtual G4double GetXmax(G4double y) override;
0074 
0075   G4double getVar(G4int);
0076   // test routine
0077   void getCoeff();
0078 
0079   G4PolarizedAnnihilationXS& operator=(const G4PolarizedAnnihilationXS& right) =
0080     delete;
0081   G4PolarizedAnnihilationXS(const G4PolarizedAnnihilationXS&) = delete;
0082 
0083  private:
0084   void TotalXS();
0085   void DefineCoefficients(const G4StokesVector& pol0,
0086                           const G4StokesVector& pol1);
0087 
0088   static constexpr G4double re2 =
0089     CLHEP::classic_electr_radius * CLHEP::classic_electr_radius;
0090 
0091   // - part depending on the polarization of the final positron
0092   G4ThreeVector fPhi2;
0093 
0094   // - part depending on the polarization of the final electron
0095   G4ThreeVector fPhi3;
0096 
0097   G4double polxx, polyy, polzz, polxz, polzx, polxy, polyx, polyz, polzy;
0098 
0099   // - unpolarised + part depending on the polarization of the initial pair
0100   G4double fPhi0;
0101 
0102   G4double fDice;
0103   G4double fPolXS, fUnpXS;
0104   G4double ISPxx, ISPyy, ISPzz, ISPnd;
0105 };
0106 #endif