source: trunk/source/processes/hadronic/cross_sections/include/G4CrossSectionDataStore.hh@ 1292

Last change on this file since 1292 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 3.8 KB
Line 
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// $Id: G4CrossSectionDataStore.hh,v 1.14 2009/01/24 11:54:47 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4CrossSectionDataStore
35//
36//
37// Modifications:
38// 23.01.2009 V.Ivanchenko move constructor and destructor to source,
39// use STL vector instead of C-array
40//
41
42// Class Description
43// This is the class to which cross section data sets may be registered.
44// An instance of it is contained in each hadronic process, allowing
45// the use of the AddDataSet() method to tailor the cross sections to
46// your application.
47// Class Description - End
48
49#ifndef G4CrossSectionDataStore_h
50#define G4CrossSectionDataStore_h 1
51
52#include "G4ParticleDefinition.hh"
53#include "G4DynamicParticle.hh"
54#include "G4Isotope.hh"
55#include "G4Element.hh"
56#include "G4Material.hh"
57#include "G4VCrossSectionDataSet.hh"
58#include <vector>
59
60class G4Nucleus;
61
62class G4CrossSectionDataStore
63{
64public:
65
66 G4CrossSectionDataStore();
67
68 ~G4CrossSectionDataStore();
69
70 G4double GetCrossSection(const G4DynamicParticle*,
71 const G4Element*, G4double aTemperature);
72
73 G4double GetCrossSection(const G4DynamicParticle*,
74 const G4Isotope*, G4double aTemperature);
75
76 G4double GetCrossSection(const G4DynamicParticle*,
77 G4double Z, G4double A, G4double aTemperature);
78
79 // to replace GetMicroscopicCrossSection
80 G4double GetCrossSection(const G4DynamicParticle*, const G4Material*);
81
82 //std::pair<G4double/*Z*/, G4double/*A*/>
83 // SelectRandomIsotope(const G4DynamicParticle*, const G4Material*);
84
85 G4Element* SampleZandA(const G4DynamicParticle*, const G4Material*,
86 G4Nucleus& target);
87
88 void AddDataSet(G4VCrossSectionDataSet*);
89
90 void BuildPhysicsTable(const G4ParticleDefinition&);
91
92 void DumpPhysicsTable(const G4ParticleDefinition&);
93
94 inline void SetVerboseLevel(G4int value)
95 {
96 verboseLevel = value;
97 }
98
99 inline G4int GetVerboseLevel()
100 {
101 return verboseLevel;
102 }
103
104private:
105
106 G4VCrossSectionDataSet* whichDataSetInCharge(const G4DynamicParticle*,
107 const G4Element*);
108
109 std::vector<G4VCrossSectionDataSet*> DataSetList;
110 G4int NDataSetList;
111 G4int verboseLevel;
112};
113
114#endif
Note: See TracBrowser for help on using the repository browser.