Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VXRayModel.cc
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 file
29//
30// File name: G4VXRayModel
31//
32// Author: Vladimir Ivanchenko
33//
34// Creation date: 28.04.2025
35//
36// -------------------------------------------------------------------
37//
38
39#include "G4VXRayModel.hh"
40#include "G4LogicalVolume.hh"
41#include "G4LossTableManager.hh"
43#include "G4DynamicParticle.hh"
45#include "G4Track.hh"
46#include <iostream>
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
52 : pName(nam)
53{
54 Register();
55}
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
60 : pLogicalVolumes(nullptr),
61 pCurrentLV(nullptr),
62 pCurrentTrack(nullptr),
63 pBetaMin(1.0),
64 pPreStepBeta(0.0),
65 pMaxBetaChange(0.1),
66 pMaxPhotons(100),
67 pTypeInt(0),
68 pVerbose(1),
69 isMaster(false),
70 pName(right.pName)
71{
72 Register();
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
78{
79 pEmManager->DeRegister(this);
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
84void G4VXRayModel::Register()
85{
87 pEmManager->Register(this);
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
93G4double G4VXRayModel::Initialise(std::vector<const G4LogicalVolume*>* ptr)
94{
95 // definition of logical volumes
96 pLogicalVolumes = ptr;
97 if (nullptr != ptr) { nVolumes = (G4int)(ptr->size()); }
98 auto params = G4OpticalParameters::Instance();
99
100 // these parameters used by Cerenkov models
101 pMaxBetaChange = params->GetCerenkovMaxBetaChange();
102 pMaxPhotons = params->GetCerenkovMaxPhotonsPerStep();
103
104 // common initialisation
105 pVerbose = params->GetCerenkovVerboseLevel();
107
108 // needed for Cerenkov process
109 return pBetaMin;
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
119G4bool G4VXRayModel::StepLimit(std::size_t idx, const G4Track& track,
120 G4double preStepBeta, G4double& limit)
121{
122 if (preStepBeta <= pBetaMin) { return false; }
123 pCurrentLV = nullptr;
124 if (idx < nVolumes) {
125 pIndex = idx;
126 pCurrentLV = (*pLogicalVolumes)[idx];
127 }
128 pCurrentTrack = &track;
129 pPreStepBeta = preStepBeta;
130 return StepLimitForVolume(limit);
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
142void G4VXRayModel::SampleXRays(std::vector<G4Track*>&, const G4Step&)
143{}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
147void G4VXRayModel::ModelDescription(std::ostream& outFile) const
148{
149 outFile << "The description for this model has not been written yet.\n";
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
static G4OpticalParameters * Instance()
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()
virtual void InitialiseModel()
virtual void ModelDescription(std::ostream &outFile) const
std::size_t pIndex
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