Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CerenkovQuasiTrackInfo Class Reference

#include <G4CerenkovQuasiTrackInfo.hh>

Inheritance diagram for G4CerenkovQuasiTrackInfo:

Public Member Functions

 G4CerenkovQuasiTrackInfo (const G4QuasiOpticalData &data, G4double pre_num_photons, G4double post_num_photons)
 ~G4CerenkovQuasiTrackInfo () override=default
void * operator new (size_t)
void operator delete (void *aCerenkovATI)
 G4CerenkovQuasiTrackInfo (const G4CerenkovQuasiTrackInfo &)=default
G4CerenkovQuasiTrackInfooperator= (const G4CerenkovQuasiTrackInfo &)=default
void Print () const override
G4QuasiOpticalData GetQuasiOpticalData () const
G4double GetPreNumPhotons () const
G4double GetPostNumPhotons () const
Public Member Functions inherited from G4VAuxiliaryTrackInformation
 G4VAuxiliaryTrackInformation ()
virtual ~G4VAuxiliaryTrackInformation ()

Static Public Member Functions

static G4CerenkovQuasiTrackInfoCast (const G4VAuxiliaryTrackInformation *const)

Detailed Description

Definition at line 44 of file G4CerenkovQuasiTrackInfo.hh.

Constructor & Destructor Documentation

◆ G4CerenkovQuasiTrackInfo() [1/2]

G4CerenkovQuasiTrackInfo::G4CerenkovQuasiTrackInfo ( const G4QuasiOpticalData & data,
G4double pre_num_photons,
G4double post_num_photons )
explicit

Definition at line 34 of file G4CerenkovQuasiTrackInfo.cc.

38 , fQuasiOpticalData(aOpticalData)
39 , fPreNumPhotons(aPreNumPhotons)
40 , fPostNumPhotons(aPostNumPhotons)
41{}

Referenced by Cast(), G4CerenkovQuasiTrackInfo(), operator delete(), and operator=().

◆ ~G4CerenkovQuasiTrackInfo()

G4CerenkovQuasiTrackInfo::~G4CerenkovQuasiTrackInfo ( )
overridedefault

◆ G4CerenkovQuasiTrackInfo() [2/2]

G4CerenkovQuasiTrackInfo::G4CerenkovQuasiTrackInfo ( const G4CerenkovQuasiTrackInfo & )
default

Member Function Documentation

◆ Cast()

G4CerenkovQuasiTrackInfo * G4CerenkovQuasiTrackInfo::Cast ( const G4VAuxiliaryTrackInformation * const aATI)
static

Definition at line 48 of file G4CerenkovQuasiTrackInfo.cc.

50{
51 G4CerenkovQuasiTrackInfo* CATI = nullptr;
52 if(aATI != nullptr)
53 {
54 // No change will be done to the pointer and to the pointed data
55 auto temp = const_cast<G4VAuxiliaryTrackInformation*>(aATI);
56 CATI = dynamic_cast<G4CerenkovQuasiTrackInfo*>(temp);
57 }
58 return CATI;
59}
G4CerenkovQuasiTrackInfo(const G4QuasiOpticalData &data, G4double pre_num_photons, G4double post_num_photons)

◆ GetPostNumPhotons()

G4double G4CerenkovQuasiTrackInfo::GetPostNumPhotons ( ) const
inline

Definition at line 65 of file G4CerenkovQuasiTrackInfo.hh.

65{ return fPostNumPhotons; }

◆ GetPreNumPhotons()

G4double G4CerenkovQuasiTrackInfo::GetPreNumPhotons ( ) const
inline

Definition at line 64 of file G4CerenkovQuasiTrackInfo.hh.

64{ return fPreNumPhotons; }

◆ GetQuasiOpticalData()

G4QuasiOpticalData G4CerenkovQuasiTrackInfo::GetQuasiOpticalData ( ) const
inline

Definition at line 63 of file G4CerenkovQuasiTrackInfo.hh.

63{ return fQuasiOpticalData; }

◆ operator delete()

void G4CerenkovQuasiTrackInfo::operator delete ( void * aCerenkovATI)
inline

Definition at line 101 of file G4CerenkovQuasiTrackInfo.hh.

102{
103 aCerenkovATIAllocator()->FreeSingle((G4CerenkovQuasiTrackInfo*) aCerenkovATI);
104}
G4DLLIMPORT G4Allocator< G4CerenkovQuasiTrackInfo > *& aCerenkovATIAllocator()

◆ operator new()

void * G4CerenkovQuasiTrackInfo::operator new ( size_t )
inline

Definition at line 92 of file G4CerenkovQuasiTrackInfo.hh.

93{
94 if(aCerenkovATIAllocator() == nullptr)
95 {
96 aCerenkovATIAllocator() = new G4Allocator<G4CerenkovQuasiTrackInfo>;
97 }
98 return (void*) aCerenkovATIAllocator()->MallocSingle();
99}

◆ operator=()

G4CerenkovQuasiTrackInfo & G4CerenkovQuasiTrackInfo::operator= ( const G4CerenkovQuasiTrackInfo & )
default

◆ Print()

void G4CerenkovQuasiTrackInfo::Print ( ) const
overridevirtual

Reimplemented from G4VAuxiliaryTrackInformation.

Definition at line 43 of file G4CerenkovQuasiTrackInfo.cc.

44{
45 G4cout << "Auxiliary track information for a Cerenkov step" << G4endl;
46}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

The documentation for this class was generated from the following files: