//////////////////////////////////////////////////////////////////// /// \class RAT::VertexGen_Flasher /// /// \brief Vertex generator for flashers /// /// \author John Walker -- contact person /// /// REVISION HISTORY /// 2015-05-28 : J. Walker - New file /// 2016-03-15 : J. Walker - Functionality to replicate SNO flasher distribution. /// 2016-05-04 : J. Walker - Set upper limit on number of photons. /// 2016-05-04 : J. Walker - Added 350nm photon distribution. /// 2017-02-10 : N. Barros - Added BeginOfRun /// /// \details This generator simulates photon emissions from PMTs. User /// RatDB settings for intensity, angular, timing and wavelength distributions. //////////////////////////////////////////////////////////////////// #ifndef __RAT_VertexGen_Flasher__ #define __RAT_VertexGen_Flasher__ #include #include #include #include #include class G4ParticleDefinition; namespace RAT { class VertexGen_Flasher : public GLG4VertexGen { public: VertexGen_Flasher( const char *arg_dbname="vertexflasher" ); virtual ~VertexGen_Flasher(){ }; void Init(); virtual void BeginOfRun(); virtual void GeneratePrimaryVertex( G4Event *event, G4ThreeVector &dx, G4double dt ); virtual void SetDirection( G4ThreeVector &dir ); G4double SampleRandom( std::vector myx, std::vector myy, std::vector mycum ); /** interface for getting and setting the state*/ virtual void SetState( G4String state ); virtual G4String GetState(); private: /** Only call the initialisation once - if this flag isn't set */ G4bool fInitDB; /** Flasher direction */ G4double fDirX, fDirY, fDirZ; /** Distribution type */ G4String fFlasherDist; /** For "uniform" distribution: range of photons per Flasher pulse */ G4int fMinPhotonsPerEvent; G4int fMaxPhotonsPerEvent; /** For "uniform" distribution: \n true = number of photons per pulse is selected from poisson distribution with mean=photons_per_event, \n false = number of photons per pulse=photons_per_event */ G4bool fPoisson; /** For "SNO" distribution: */ G4double fNphotonNhitParam0; G4double fNphotonNhitParam1; G4double fNphotonNhitParam2; G4double fSNODistParam0; G4double fSNODistParam1; G4double fSNODistParam2; G4double fSNODistParam3; G4double fSNODistParam4; G4double fSNODistParam5; G4double fSNODistParam6; G4double fSNODistParam7; G4double fSNODistParam8; G4double fSNODistParam9; G4double fSNODistBound0; G4double fSNODistBound1; G4double fSNODistBound2; G4double fSNODistBound3; G4int fSNOMaxNPhoton; /** Wavelength mode for flashers */ /** flag set to single wavelength value */ G4double fMono_wl; /** flag set to temperature for blackbody radiation spectrum */ G4int fBlackbody_T; /** Flasher angular distribution - string tells you what distribution to read */ G4String fFlasherWavelength; /** Time mode for flashers (1 for "mono", 0 for the distribution) */ G4int fMono_time; /** Flasher time distribution - string tells you what distribution to read */ G4String fFlasherTime; /** Angle mode for flashers */ /** flag set to 1 for isotropic dispersion of the beam */ G4int fIso_angle; /** flag set to 1 if there is to be no dispersion of the beam */ G4int fZero_angle; /** Flasher angular distribution - string tells you what distribution to read */ G4String fFlasherAngle; /** Vector for specified angles */ std::vector fang; /** Vector for cumulative angular intensity distribution */ std::vector fangCumMag; /** Vector for angular intensity */ std::vector fangInt; /** Vector for specified times */ std::vector ftime; /** Vector for cumulative time intensity distribution */ std::vector ftimeCumMag; /** Vector for time intensity */ std::vector ftimeInt; /** Vector for specified wavelengths */ std::vector fwave; /** Vector for cumulative wave intensity distribution */ std::vector fwaveCumMag; /** Vector for wave intensity */ std::vector fwaveInt; /** Vector for wave probability */ std::vector fwaveProb; /** Vector for wave cumulative probability */ std::vector fwaveCumProb; /** Limits in angular distribution */ G4double fMinAngle, fMaxAngle; /** Limits on wavelength distribution */ G4double fMinWl, fMaxWl; /** Limits on time distribution */ G4double fMinTime,fMaxTime; /** Factor for converting angles (random generator uses radians so if distribution set in degrees = 360/twopi, otherwise =1*/ G4double fConvertAngle; /** Particles must be defined as optical photons */ G4ParticleDefinition *fPhotonDef; }; } // namespace RAT #endif // __RAT_VertexGen_Flasher__