Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/include/G4CrossSectionDataSet.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 // Author: Riccardo Capra <capra@ge.infn.it>
 29 //
 30 // History:
 31 // -----------
 32 // 30 Jun 2005  RC                  Created
 33 // 14 Oct 2007  MGP                 Removed inheritance from concrete class G4ShellEMDataSet
 34 // 15 Jul 2009  N.A.Karakatsanis    New methods added for loading logarithmic data
 35 //                                  to enhance computing performance of interpolation
 36 //
 37 // -------------------------------------------------------------------
 38 
 39 // Class description:
 40 // Low Energy Electromagnetic Physics
 41 // Data set for an electromagnetic physics process
 42 // A strategy pattern is used to encapsulate algorithms for data interpolation
 43 // -------------------------------------------------------------------
 44 
 45 #ifndef  G4CrossSectionDataSet_HH
 46 #define  G4CrossSectionDataSet_HH 1
 47 
 48 #include <CLHEP/Units/SystemOfUnits.h>
 49 
 50 #include "G4ShellEMDataSet.hh"
 51 
 52 class G4CrossSectionDataSet : public G4VEMDataSet
 53 { 
 54 
 55 public:
 56   explicit G4CrossSectionDataSet(G4VDataSetAlgorithm* algo, 
 57          G4double xUnit=CLHEP::MeV, 
 58          G4double dataUnit=CLHEP::barn);
 59 
 60   virtual ~G4CrossSectionDataSet();
 61 
 62   virtual G4double FindValue(G4double e, G4int componentId=0) const;
 63   
 64   virtual void PrintData(void) const;
 65 
 66   virtual const G4VEMDataSet*  GetComponent(G4int componentId) const 
 67   { return components[componentId]; }
 68 
 69   virtual void AddComponent(G4VEMDataSet* dataSet) 
 70   { components.push_back(dataSet); }
 71 
 72   virtual size_t NumberOfComponents(void) const 
 73   { return components.size(); }
 74 
 75   virtual const G4DataVector& GetEnergies(G4int componentId) const 
 76   { return GetComponent(componentId)->GetEnergies(0); }
 77 
 78   virtual const G4DataVector& GetData(G4int componentId) const 
 79   { return GetComponent(componentId)->GetData(0); }
 80 
 81   virtual const G4DataVector& GetLogEnergies(G4int componentId) const 
 82   { return GetComponent(componentId)->GetLogEnergies(0); }
 83 
 84   virtual const G4DataVector& GetLogData(G4int componentId) const 
 85   { return GetComponent(componentId)->GetLogData(0); }
 86 
 87   virtual void SetEnergiesData(G4DataVector* x, G4DataVector* values, G4int componentId);
 88 
 89   virtual void SetLogEnergiesData(G4DataVector* x,
 90                                   G4DataVector* values,
 91                                   G4DataVector* log_x, 
 92                                   G4DataVector* log_values,
 93                                   G4int componentId);
 94 
 95   virtual G4bool LoadData(const G4String & argFileName);
 96   virtual G4bool LoadNonLogData(const G4String & argFileName);
 97 
 98   virtual G4bool SaveData(const G4String & argFileName) const;
 99  
100   virtual G4double RandomSelect(G4int /*componentId */) const { return -1.; };
101 
102 private:
103 
104   G4String FullFileName(const G4String & argFileName) const;
105 
106   // Hide copy constructor and assignment operator 
107   explicit G4CrossSectionDataSet();
108   G4CrossSectionDataSet(const G4CrossSectionDataSet & copy) = delete;
109   G4CrossSectionDataSet& operator=(const G4CrossSectionDataSet & right) = delete;
110 
111   G4double GetUnitEnergies() const { return unitEnergies; }
112   G4double GetUnitData() const { return unitData; }
113   const G4VDataSetAlgorithm* GetAlgorithm() const { return algorithm; }
114    
115   void CleanUpComponents(void);
116 
117   std::vector<G4VEMDataSet*> components;          // Owned pointers
118 
119   G4VDataSetAlgorithm* algorithm;           // Owned pointer 
120   
121   G4double unitEnergies;
122   G4double unitData; 
123 
124   G4int z = 0;
125 };
126 #endif /* G4CrossSectionDataSet_HH */
127