Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VXRayModel.hh
Go to the documentation of this file.
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//
31// File name: G4VXRayModel
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 28.04.2025
36//
37//
38// Class Description:
39//
40// Abstract interface to a X-Ray production model
41
42// -------------------------------------------------------------------
43//
44
45#ifndef G4VXRayModel_h
46#define G4VXRayModel_h 1
47
48#include "globals.hh"
49#include <vector>
50
51class G4LogicalVolume;
55class G4Track;
56class G4Step;
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
61{
62
63public:
64
65 explicit G4VXRayModel(const G4String& nam);
66
68
69 virtual ~G4VXRayModel();
70
71 // return minimal beta for Cerenkov in these volumes
72 // the method may be used for other X-Ray models
73 G4double Initialise(std::vector<const G4LogicalVolume*>*);
74 virtual void InitialiseModel();
75
76 // check applicability and propose step limit, which may be DBL_MAX
77 // if these methods return "false", then sampling of X-rays not possible
78 G4bool StepLimit(std::size_t idx, const G4Track&,
79 G4double preStepBeta, G4double& limit);
80 virtual G4bool StepLimitForVolume(G4double& limit);
81
82 // sampling is called if StepLimit(..) returns "true"
83 // produced X-rays are inside vector out
84 // each photon has time, position, and other parameters within the step
85 virtual void SampleXRays(std::vector<G4Track*>& out, const G4Step&);
86
87 // for automatic documentation
88 virtual void ModelDescription(std::ostream& outFile) const;
89
90 const G4String& GetName() const { return pName; };
91
92 G4int GetType() const { return pTypeInt; };
93
94 // hide assignment operator
95 G4VXRayModel& operator=(const G4VXRayModel& right) = delete;
96
97private:
98
99 void Register();
100
101protected:
102
103 std::vector<const G4LogicalVolume*>* pLogicalVolumes{nullptr};
106 const G4Track* pCurrentTrack{nullptr};
110
112 std::size_t pIndex{0};
113 std::size_t nVolumes{0};
114
119};
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
122
123#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
const G4String & GetName() const
std::vector< const G4LogicalVolume * > * pLogicalVolumes
const G4LogicalVolume * pCurrentLV
virtual void SampleXRays(std::vector< G4Track * > &out, const G4Step &)
G4double pMaxBetaChange
G4VXRayModel(const G4String &nam)
virtual ~G4VXRayModel()
G4int GetType() const
virtual void InitialiseModel()
virtual void ModelDescription(std::ostream &outFile) const
std::size_t pIndex
G4VXRayModel & operator=(const G4VXRayModel &right)=delete
G4bool StepLimit(std::size_t idx, const G4Track &, G4double preStepBeta, G4double &limit)
G4LossTableManager * pEmManager
G4double Initialise(std::vector< const G4LogicalVolume * > *)
G4double pBetaMin
virtual G4bool StepLimitForVolume(G4double &limit)
const G4Track * pCurrentTrack
std::size_t nVolumes
G4double pPreStepBeta
const G4String pName