Geant4  10.02.p03
G4Scintillation Class Reference

#include <G4Scintillation.hh>

Inheritance diagram for G4Scintillation:
Collaboration diagram for G4Scintillation:

Public Member Functions

 G4Scintillation (const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
 
 ~G4Scintillation ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4VParticleChange * AtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4double GetScintillationYieldByParticleType (const G4Track &aTrack, const G4Step &aStep)
 
void SetTrackSecondariesFirst (const G4bool state)
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
G4bool GetFiniteRiseTime () const
 
void SetScintillationYieldFactor (const G4double yieldfactor)
 
G4double GetScintillationYieldFactor () const
 
void SetScintillationExcitationRatio (const G4double ratio)
 
G4double GetScintillationExcitationRatio () const
 
G4PhysicsTableGetFastIntegralTable () const
 
G4PhysicsTableGetSlowIntegralTable () const
 
void AddSaturation (G4EmSaturation *)
 
void RemoveSaturation ()
 
G4EmSaturationGetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void DumpPhysicsTable () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

void BuildThePhysicsTable ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTablefFastIntegralTable
 
G4PhysicsTablefSlowIntegralTable
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Private Member Functions

 G4Scintillation (const G4Scintillation &right)
 
G4Scintillationoperator= (const G4Scintillation &right)
 
G4double single_exp (G4double t, G4double tau2)
 
G4double bi_exp (G4double t, G4double tau1, G4double tau2)
 
G4double sample_time (G4double tau1, G4double tau2)
 

Private Attributes

G4bool fTrackSecondariesFirst
 
G4bool fFiniteRiseTime
 
G4double fYieldFactor
 
G4double fExcitationRatio
 
G4bool fScintillationByParticleType
 
G4EmSaturationfEmSaturation
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 88 of file G4Scintillation.hh.

Constructor & Destructor Documentation

◆ G4Scintillation() [1/2]

G4Scintillation::G4Scintillation ( const G4String processName = "Scintillation",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 108 of file G4Scintillation.cc.

110  : G4VRestDiscreteProcess(processName, type) ,
111  fTrackSecondariesFirst(false),
112  fFiniteRiseTime(false),
113  fYieldFactor(1.0),
114  fExcitationRatio(1.0),
116  fEmSaturation(nullptr)
117 {
119 
120 #ifdef G4DEBUG_SCINTILLATION
121  ScintTrackEDep = 0.;
122  ScintTrackYield = 0.;
123 #endif
124 
125  fFastIntegralTable = nullptr;
126  fSlowIntegralTable = nullptr;
127 
128  if (verboseLevel>0) {
129  G4cout << GetProcessName() << " is created " << G4endl;
130  }
131 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4bool fTrackSecondariesFirst
G4PhysicsTable * fFastIntegralTable
G4EmSaturation * fEmSaturation
G4PhysicsTable * fSlowIntegralTable
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4double fExcitationRatio
#define G4endl
Definition: G4ios.hh:61
G4bool fScintillationByParticleType
Here is the call graph for this function:

◆ ~G4Scintillation()

G4Scintillation::~G4Scintillation ( )

Definition at line 137 of file G4Scintillation.cc.

138 {
139  if (fFastIntegralTable != nullptr) {
141  delete fFastIntegralTable;
142  }
143  if (fSlowIntegralTable != nullptr) {
145  delete fSlowIntegralTable;
146  }
147 }
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
void clearAndDestroy()
Here is the call graph for this function:

◆ G4Scintillation() [2/2]

G4Scintillation::G4Scintillation ( const G4Scintillation right)
private

Member Function Documentation

◆ AddSaturation()

void G4Scintillation::AddSaturation ( G4EmSaturation sat)

Definition at line 674 of file G4Scintillation.cc.

675 {
676  fEmSaturation = sat;
677 }
G4EmSaturation * fEmSaturation
Here is the caller graph for this function:

◆ AtRestDoIt()

G4VParticleChange * G4Scintillation::AtRestDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 173 of file G4Scintillation.cc.

178 {
179  return G4Scintillation::PostStepDoIt(aTrack, aStep);
180 }
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Here is the call graph for this function:

◆ bi_exp()

G4double G4Scintillation::bi_exp ( G4double  t,
G4double  tau1,
G4double  tau2 
)
inlineprivate

Definition at line 333 of file G4Scintillation.hh.

334 {
335  return std::exp(-1.0*t/tau2)*(1-std::exp(-1.0*t/tau1))/tau2/tau2*(tau1+tau2);
336 }
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4Scintillation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 153 of file G4Scintillation.cc.

154 {
155  if (fFastIntegralTable != nullptr) {
157  delete fFastIntegralTable;
158  fFastIntegralTable = nullptr;
159  }
160  if (fSlowIntegralTable != nullptr) {
162  delete fSlowIntegralTable;
163  fSlowIntegralTable = nullptr;
164  }
166 }
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
void clearAndDestroy()
Here is the call graph for this function:

◆ BuildThePhysicsTable()

void G4Scintillation::BuildThePhysicsTable ( )
protected

Definition at line 489 of file G4Scintillation.cc.

490 {
491  if (fFastIntegralTable && fSlowIntegralTable) return;
492 
493  const G4MaterialTable* theMaterialTable =
495  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
496 
497  // create new physics table
498 
500  new G4PhysicsTable(numOfMaterials);
502  new G4PhysicsTable(numOfMaterials);
503 
504  // loop for materials
505 
506  for (G4int i=0 ; i < numOfMaterials; i++)
507  {
508  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
510  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
512 
513  // Retrieve vector of scintillation wavelength intensity for
514  // the material from the material's optical properties table.
515 
516  G4Material* aMaterial = (*theMaterialTable)[i];
517 
518  G4MaterialPropertiesTable* aMaterialPropertiesTable =
519  aMaterial->GetMaterialPropertiesTable();
520 
521  if (aMaterialPropertiesTable) {
522 
523  G4MaterialPropertyVector* theFastLightVector =
524  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
525 
526  if (theFastLightVector) {
527 
528  // Retrieve the first intensity point in vector
529  // of (photon energy, intensity) pairs
530 
531  G4double currentIN = (*theFastLightVector)[0];
532 
533  if (currentIN >= 0.0) {
534 
535  // Create first (photon energy, Scintillation
536  // Integral pair
537 
538  G4double currentPM = theFastLightVector->Energy(0);
539 
540  G4double currentCII = 0.0;
541 
542  aPhysicsOrderedFreeVector->
543  InsertValues(currentPM , currentCII);
544 
545  // Set previous values to current ones prior to loop
546 
547  G4double prevPM = currentPM;
548  G4double prevCII = currentCII;
549  G4double prevIN = currentIN;
550 
551  // loop over all (photon energy, intensity)
552  // pairs stored for this material
553 
554  for (size_t ii = 1;
555  ii < theFastLightVector->GetVectorLength();
556  ++ii)
557  {
558  currentPM = theFastLightVector->Energy(ii);
559  currentIN = (*theFastLightVector)[ii];
560 
561  currentCII = 0.5 * (prevIN + currentIN);
562 
563  currentCII = prevCII +
564  (currentPM - prevPM) * currentCII;
565 
566  aPhysicsOrderedFreeVector->
567  InsertValues(currentPM, currentCII);
568 
569  prevPM = currentPM;
570  prevCII = currentCII;
571  prevIN = currentIN;
572  }
573 
574  }
575  }
576 
577  G4MaterialPropertyVector* theSlowLightVector =
578  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
579 
580  if (theSlowLightVector) {
581 
582  // Retrieve the first intensity point in vector
583  // of (photon energy, intensity) pairs
584 
585  G4double currentIN = (*theSlowLightVector)[0];
586 
587  if (currentIN >= 0.0) {
588 
589  // Create first (photon energy, Scintillation
590  // Integral pair
591 
592  G4double currentPM = theSlowLightVector->Energy(0);
593 
594  G4double currentCII = 0.0;
595 
596  bPhysicsOrderedFreeVector->
597  InsertValues(currentPM , currentCII);
598 
599  // Set previous values to current ones prior to loop
600 
601  G4double prevPM = currentPM;
602  G4double prevCII = currentCII;
603  G4double prevIN = currentIN;
604 
605  // loop over all (photon energy, intensity)
606  // pairs stored for this material
607 
608  for (size_t ii = 1;
609  ii < theSlowLightVector->GetVectorLength();
610  ++ii)
611  {
612  currentPM = theSlowLightVector->Energy(ii);
613  currentIN = (*theSlowLightVector)[ii];
614 
615  currentCII = 0.5 * (prevIN + currentIN);
616 
617  currentCII = prevCII +
618  (currentPM - prevPM) * currentCII;
619 
620  bPhysicsOrderedFreeVector->
621  InsertValues(currentPM, currentCII);
622 
623  prevPM = currentPM;
624  prevCII = currentCII;
625  prevIN = currentIN;
626  }
627 
628  }
629  }
630  }
631 
632  // The scintillation integral(s) for a given material
633  // will be inserted in the table(s) according to the
634  // position of the material in the material table.
635 
636  fFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
637  fSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
638 
639  }
640 }
G4MaterialPropertyVector * GetProperty(const char *key)
G4PhysicsTable * fFastIntegralTable
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
G4PhysicsTable * fSlowIntegralTable
int G4int
Definition: G4Types.hh:78
size_t GetVectorLength() const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:596
void insertAt(size_t, G4PhysicsVector *)
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DumpPhysicsTable()

void G4Scintillation::DumpPhysicsTable ( ) const
inline

Definition at line 301 of file G4Scintillation.hh.

302 {
303  if (fFastIntegralTable) {
304  G4int PhysicsTableSize = fFastIntegralTable->entries();
306 
307  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
308  {
310  v->DumpValues();
311  }
312  }
313 
314  if (fSlowIntegralTable) {
315  G4int PhysicsTableSize = fSlowIntegralTable->entries();
317 
318  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
319  {
321  v->DumpValues();
322  }
323  }
324 }
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
int G4int
Definition: G4Types.hh:78
size_t entries() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetFastIntegralTable()

G4PhysicsTable * G4Scintillation::GetFastIntegralTable ( ) const
inline

Definition at line 295 of file G4Scintillation.hh.

296 {
297  return fFastIntegralTable;
298 }
G4PhysicsTable * fFastIntegralTable

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 271 of file G4Scintillation.hh.

272 {
273  return fFiniteRiseTime;
274 }

◆ GetMeanFreePath()

G4double G4Scintillation::GetMeanFreePath ( const G4Track &  aTrack,
G4double  ,
G4ForceCondition *  condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 688 of file G4Scintillation.cc.

691 {
692  *condition = StronglyForced;
693 
694  return DBL_MAX;
695 
696 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ GetMeanLifeTime()

G4double G4Scintillation::GetMeanLifeTime ( const G4Track &  aTrack,
G4ForceCondition *  condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 702 of file G4Scintillation.cc.

704 {
705  *condition = Forced;
706 
707  return DBL_MAX;
708 
709 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ GetSaturation()

G4EmSaturation* G4Scintillation::GetSaturation ( ) const
inline

Definition at line 198 of file G4Scintillation.hh.

198 { return fEmSaturation; }
G4EmSaturation * fEmSaturation
Here is the call graph for this function:

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 205 of file G4Scintillation.hh.

206  { return fScintillationByParticleType; }
G4bool fScintillationByParticleType
Here is the call graph for this function:

◆ GetScintillationExcitationRatio()

G4double G4Scintillation::GetScintillationExcitationRatio ( ) const
inline

Definition at line 283 of file G4Scintillation.hh.

284 {
285  return fExcitationRatio;
286 }
G4double fExcitationRatio

◆ GetScintillationYieldByParticleType()

G4double G4Scintillation::GetScintillationYieldByParticleType ( const G4Track &  aTrack,
const G4Step &  aStep 
)

Definition at line 734 of file G4Scintillation.cc.

735 {
737  // Get the scintillation yield vector //
739 
740  G4ParticleDefinition *pDef = aTrack.GetDynamicParticle()->GetDefinition();
741 
742  G4MaterialPropertyVector *Scint_Yield_Vector = nullptr;
743 
744  G4MaterialPropertiesTable *aMaterialPropertiesTable
745  = aTrack.GetMaterial()->GetMaterialPropertiesTable();
746 
747  // Get the G4MaterialPropertyVector containing the scintillation
748  // yield as a function of the energy deposited and particle type
749 
750  // Protons
751  if(pDef==G4Proton::ProtonDefinition())
752  Scint_Yield_Vector = aMaterialPropertiesTable->
753  GetProperty("PROTONSCINTILLATIONYIELD");
754 
755  // Deuterons
756  else if(pDef==G4Deuteron::DeuteronDefinition())
757  Scint_Yield_Vector = aMaterialPropertiesTable->
758  GetProperty("DEUTERONSCINTILLATIONYIELD");
759 
760  // Tritons
761  else if(pDef==G4Triton::TritonDefinition())
762  Scint_Yield_Vector = aMaterialPropertiesTable->
763  GetProperty("TRITONSCINTILLATIONYIELD");
764 
765  // Alphas
766  else if(pDef==G4Alpha::AlphaDefinition())
767  Scint_Yield_Vector = aMaterialPropertiesTable->
768  GetProperty("ALPHASCINTILLATIONYIELD");
769 
770  // Ions (particles derived from G4VIon and G4Ions) and recoil ions
771  // below the production cut from neutrons after hElastic
772  else if(pDef->GetParticleType()== "nucleus" ||
774  Scint_Yield_Vector = aMaterialPropertiesTable->
775  GetProperty("IONSCINTILLATIONYIELD");
776 
777  // Electrons (must also account for shell-binding energy
778  // attributed to gamma from standard photoelectric effect)
779  else if(pDef==G4Electron::ElectronDefinition() ||
780  pDef==G4Gamma::GammaDefinition())
781  Scint_Yield_Vector = aMaterialPropertiesTable->
782  GetProperty("ELECTRONSCINTILLATIONYIELD");
783 
784  // Default for particles not enumerated/listed above
785  else
786  Scint_Yield_Vector = aMaterialPropertiesTable->
787  GetProperty("ELECTRONSCINTILLATIONYIELD");
788 
789  // If the user has specified none of the above particles then the
790  // default is the electron scintillation yield
791  if(!Scint_Yield_Vector)
792  Scint_Yield_Vector = aMaterialPropertiesTable->
793  GetProperty("ELECTRONSCINTILLATIONYIELD");
794 
795  // Throw an exception if no scintillation yield vector is found
796  if (!Scint_Yield_Vector) {
798  ed << "\nG4Scintillation::PostStepDoIt(): "
799  << "Request for scintillation yield for energy deposit and particle\n"
800  << "type without correct entry in MaterialPropertiesTable.\n"
801  << "ScintillationByParticleType requires at minimum that \n"
802  << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
803  << G4endl;
804  G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable";
805  G4Exception("G4Scintillation::PostStepDoIt","Scint01",
806  FatalException,ed,comments);
807  }
808 
810  // Calculate the scintillation light //
812  // To account for potential nonlinearity and scintillation photon
813  // density along the track, light (L) is produced according to:
814  //
815  // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
816 
817  G4double ScintillationYield = 0.;
818 
819  G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
820  G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
821 
822  if(PreStepKineticEnergy <= Scint_Yield_Vector->GetMaxEnergy()){
823  G4double Yield1 = Scint_Yield_Vector->Value(PreStepKineticEnergy);
824  G4double Yield2 = Scint_Yield_Vector->
825  Value(PreStepKineticEnergy - StepEnergyDeposit);
826  ScintillationYield = Yield1 - Yield2;
827  } else {
829  ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
830  << "for scintillation light yield above the available energy range\n"
831  << "specifed in G4MaterialPropertiesTable. A linear interpolation\n"
832  << "will be performed to compute the scintillation light yield using\n"
833  << "(L_max / E_max) as the photon yield per unit energy."
834  << G4endl;
835  G4String cmt = "\nScintillation yield may be unphysical!\n";
836  G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
837  "Scint03", JustWarning, ed, cmt);
838 
839  G4double LinearYield = Scint_Yield_Vector->GetMaxValue()
840  / Scint_Yield_Vector->GetMaxEnergy();
841 
842  // Units: [# scintillation photons]
843  ScintillationYield = LinearYield * StepEnergyDeposit;
844  }
845 
846 #ifdef G4DEBUG_SCINTILLATION
847 
848  // Increment track aggregators
849  ScintTrackYield += ScintillationYield;
850  ScintTrackEDep += StepEnergyDeposit;
851 
852  G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
853  << "--\n"
854  << "-- Name = " << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
855  << "-- TrackID = " << aTrack.GetTrackID() << "\n"
856  << "-- ParentID = " << aTrack.GetParentID() << "\n"
857  << "-- Current KE = " << aTrack.GetKineticEnergy()/MeV << " MeV\n"
858  << "-- Step EDep = " << aStep.GetTotalEnergyDeposit()/MeV << " MeV\n"
859  << "-- Track EDep = " << ScintTrackEDep/MeV << " MeV\n"
860  << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy()/MeV << " MeV\n"
861  << "-- Step yield = " << ScintillationYield << " photons\n"
862  << "-- Track yield = " << ScintTrackYield << " photons\n"
863  << G4endl;
864 
865  // The track has terminated within or has left the scintillator volume
866  if( (aTrack.GetTrackStatus() == fStopButAlive) or
867  (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) ){
868 
869  // Reset aggregators for the next track
870  ScintTrackEDep = 0.;
871  ScintTrackYield = 0.;
872  }
873 
874 #endif
875 
876  return ScintillationYield;
877 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:90
static const double MeV
Definition: G4SIunits.hh:211
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
const G4String & GetParticleType() const
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetMaxEnergy() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetScintillationYieldFactor()

G4double G4Scintillation::GetScintillationYieldFactor ( ) const
inline

Definition at line 277 of file G4Scintillation.hh.

278 {
279  return fYieldFactor;
280 }

◆ GetSlowIntegralTable()

G4PhysicsTable * G4Scintillation::GetSlowIntegralTable ( ) const
inline

Definition at line 289 of file G4Scintillation.hh.

290 {
291  return fSlowIntegralTable;
292 }
G4PhysicsTable * fSlowIntegralTable

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 265 of file G4Scintillation.hh.

266 {
267  return fTrackSecondariesFirst;
268 }
G4bool fTrackSecondariesFirst

◆ IsApplicable()

G4bool G4Scintillation::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 256 of file G4Scintillation.hh.

257 {
258  if (aParticleType.GetParticleName() == "opticalphoton") return false;
259  if (aParticleType.IsShortLived()) return false;
260 
261  return true;
262 }
const G4String & GetParticleName() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4Scintillation& G4Scintillation::operator= ( const G4Scintillation right)
private

◆ PostStepDoIt()

G4VParticleChange * G4Scintillation::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 186 of file G4Scintillation.cc.

193 {
194  aParticleChange.Initialize(aTrack);
195 
196  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
197  const G4Material* aMaterial = aTrack.GetMaterial();
198 
199  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
200  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
201 
202  G4ThreeVector x0 = pPreStepPoint->GetPosition();
203  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
204  G4double t0 = pPreStepPoint->GetGlobalTime();
205 
206  G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
207 
208  G4MaterialPropertiesTable* aMaterialPropertiesTable =
209  aMaterial->GetMaterialPropertiesTable();
210  if (!aMaterialPropertiesTable)
211  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
212 
213  G4MaterialPropertyVector* Fast_Intensity =
214  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
215  G4MaterialPropertyVector* Slow_Intensity =
216  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
217 
218  if (!Fast_Intensity && !Slow_Intensity )
219  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
220 
221  G4int nscnt = 1;
222  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
223 
224  G4double ScintillationYield = 0.;
225 
226  // Scintillation depends on particle type, energy deposited
228 
229  ScintillationYield =
231 
232  // The default linear scintillation process
233  } else {
234 
235  ScintillationYield = aMaterialPropertiesTable->
236  GetConstProperty("SCINTILLATIONYIELD");
237 
238  // Units: [# scintillation photons / MeV]
239  ScintillationYield *= fYieldFactor;
240  }
241 
242  G4double ResolutionScale = aMaterialPropertiesTable->
243  GetConstProperty("RESOLUTIONSCALE");
244 
245  // Birks law saturation:
246 
247  //G4double constBirks = 0.0;
248 
249  //constBirks = aMaterial->GetIonisation()->GetBirksConstant();
250 
251  G4double MeanNumberOfPhotons;
252 
253  // Birk's correction via fEmSaturation and specifying scintillation by
254  // by particle type are physically mutually exclusive
255 
257  MeanNumberOfPhotons = ScintillationYield;
258  else if (fEmSaturation)
259  MeanNumberOfPhotons = ScintillationYield*
261  else
262  MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit;
263 
264  G4int NumPhotons;
265 
266  if (MeanNumberOfPhotons > 10.)
267  {
268  G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
269  NumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
270  }
271  else
272  {
273  NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
274  }
275 
276  if (NumPhotons <= 0)
277  {
278  // return unchanged particle and no secondaries
279 
280  aParticleChange.SetNumberOfSecondaries(0);
281 
282  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
283  }
284 
286 
287  aParticleChange.SetNumberOfSecondaries(NumPhotons);
288 
290  if (aTrack.GetTrackStatus() == fAlive )
291  aParticleChange.ProposeTrackStatus(fSuspend);
292  }
293 
295 
296  G4int materialIndex = aMaterial->GetIndex();
297 
298  // Retrieve the Scintillation Integral for this material
299  // new G4PhysicsOrderedFreeVector allocated to hold CII's
300 
301  G4int Num = NumPhotons;
302 
303  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
304 
305  G4double ScintillationTime = 0.*ns;
306  G4double ScintillationRiseTime = 0.*ns;
307  G4PhysicsOrderedFreeVector* ScintillationIntegral = nullptr;
308 
309  if (scnt == 1) {
310  if (nscnt == 1) {
311  if (Fast_Intensity) {
312  ScintillationTime = aMaterialPropertiesTable->
313  GetConstProperty("FASTTIMECONSTANT");
314  if (fFiniteRiseTime) {
315  ScintillationRiseTime = aMaterialPropertiesTable->
316  GetConstProperty("FASTSCINTILLATIONRISETIME");
317  }
318  ScintillationIntegral =
320  ((*fFastIntegralTable)(materialIndex));
321  }
322  if (Slow_Intensity) {
323  ScintillationTime = aMaterialPropertiesTable->
324  GetConstProperty("SLOWTIMECONSTANT");
325  if (fFiniteRiseTime) {
326  ScintillationRiseTime = aMaterialPropertiesTable->
327  GetConstProperty("SLOWSCINTILLATIONRISETIME");
328  }
329  ScintillationIntegral =
331  ((*fSlowIntegralTable)(materialIndex));
332  }
333  }
334  else {
335  G4double yieldRatio = aMaterialPropertiesTable->
336  GetConstProperty("YIELDRATIO");
337  if ( fExcitationRatio == 1.0 || fExcitationRatio == 0.0) {
338  Num = G4int (std::min(yieldRatio,1.0) * NumPhotons);
339  }
340  else {
341  Num = G4int (std::min(fExcitationRatio,1.0) * NumPhotons);
342  }
343  ScintillationTime = aMaterialPropertiesTable->
344  GetConstProperty("FASTTIMECONSTANT");
345  if (fFiniteRiseTime) {
346  ScintillationRiseTime = aMaterialPropertiesTable->
347  GetConstProperty("FASTSCINTILLATIONRISETIME");
348  }
349  ScintillationIntegral =
351  ((*fFastIntegralTable)(materialIndex));
352  }
353  }
354  else {
355  Num = NumPhotons - Num;
356  ScintillationTime = aMaterialPropertiesTable->
357  GetConstProperty("SLOWTIMECONSTANT");
358  if (fFiniteRiseTime) {
359  ScintillationRiseTime = aMaterialPropertiesTable->
360  GetConstProperty("SLOWSCINTILLATIONRISETIME");
361  }
362  ScintillationIntegral =
364  ((*fSlowIntegralTable)(materialIndex));
365  }
366 
367  if (!ScintillationIntegral) continue;
368 
369  // Max Scintillation Integral
370 
371  G4double CIImax = ScintillationIntegral->GetMaxValue();
372 
373  for (G4int i = 0; i < Num; i++) {
374 
375  // Determine photon energy
376 
377  G4double CIIvalue = G4UniformRand()*CIImax;
378  G4double sampledEnergy =
379  ScintillationIntegral->GetEnergy(CIIvalue);
380 
381  if (verboseLevel>1) {
382  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
383  G4cout << "CIIvalue = " << CIIvalue << G4endl;
384  }
385 
386  // Generate random photon direction
387 
388  G4double cost = 1. - 2.*G4UniformRand();
389  G4double sint = std::sqrt((1.-cost)*(1.+cost));
390 
391  G4double phi = twopi*G4UniformRand();
392  G4double sinp = std::sin(phi);
393  G4double cosp = std::cos(phi);
394 
395  G4double px = sint*cosp;
396  G4double py = sint*sinp;
397  G4double pz = cost;
398 
399  // Create photon momentum direction vector
400 
401  G4ParticleMomentum photonMomentum(px, py, pz);
402 
403  // Determine polarization of new photon
404 
405  G4double sx = cost*cosp;
406  G4double sy = cost*sinp;
407  G4double sz = -sint;
408 
409  G4ThreeVector photonPolarization(sx, sy, sz);
410 
411  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
412 
413  phi = twopi*G4UniformRand();
414  sinp = std::sin(phi);
415  cosp = std::cos(phi);
416 
417  photonPolarization = cosp * photonPolarization + sinp * perp;
418 
419  photonPolarization = photonPolarization.unit();
420 
421  // Generate a new photon:
422 
423  G4DynamicParticle* aScintillationPhoton =
425  photonMomentum);
426  aScintillationPhoton->SetPolarization
427  (photonPolarization.x(),
428  photonPolarization.y(),
429  photonPolarization.z());
430 
431  aScintillationPhoton->SetKineticEnergy(sampledEnergy);
432 
433  // Generate new G4Track object:
434 
435  G4double rand;
436 
437  if (aParticle->GetDefinition()->GetPDGCharge() != 0) {
438  rand = G4UniformRand();
439  } else {
440  rand = 1.0;
441  }
442 
443  G4double delta = rand * aStep.GetStepLength();
444  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
445  rand*(pPostStepPoint->GetVelocity()-
446  pPreStepPoint->GetVelocity())/2.);
447 
448  // emission time distribution
449  if (ScintillationRiseTime==0.0) {
450  deltaTime = deltaTime -
451  ScintillationTime * std::log( G4UniformRand() );
452  } else {
453  deltaTime = deltaTime +
454  sample_time(ScintillationRiseTime, ScintillationTime);
455  }
456 
457  G4double aSecondaryTime = t0 + deltaTime;
458 
459  G4ThreeVector aSecondaryPosition =
460  x0 + rand * aStep.GetDeltaPosition();
461 
462  G4Track* aSecondaryTrack = new G4Track(aScintillationPhoton,
463  aSecondaryTime,
464  aSecondaryPosition);
465 
466  aSecondaryTrack->SetTouchableHandle(
467  aStep.GetPreStepPoint()->GetTouchableHandle());
468  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
469 
470  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
471 
472  aParticleChange.AddSecondary(aSecondaryTrack);
473 
474  }
475  }
476 
477  if (verboseLevel>0) {
478  G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
479  << aParticleChange.GetNumberOfSecondaries() << G4endl;
480  }
481 
482  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
483 }
ThreeVector shoot(const G4int Ap, const G4int Af)
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4MaterialPropertyVector * GetProperty(const char *key)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4bool fTrackSecondariesFirst
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4PhysicsTable * fFastIntegralTable
G4EmSaturation * fEmSaturation
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4PhysicsTable * fSlowIntegralTable
size_t GetIndex() const
Definition: G4Material.hh:262
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector unit() const
static const double twopi
Definition: G4SIunits.hh:75
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void SetKineticEnergy(G4double aEnergy)
G4double GetEnergy(G4double aValue)
G4double fExcitationRatio
static G4OpticalPhoton * OpticalPhoton()
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double sample_time(G4double tau1, G4double tau2)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
#define ns
Definition: xmlparse.cc:614
G4double GetPDGCharge() const
G4double VisibleEnergyDepositionAtAStep(const G4Step *)
G4bool fScintillationByParticleType
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RemoveSaturation()

void G4Scintillation::RemoveSaturation ( )

Definition at line 679 of file G4Scintillation.cc.

680 {
681  fEmSaturation = nullptr;
682 }
G4EmSaturation * fEmSaturation
Here is the caller graph for this function:

◆ sample_time()

G4double G4Scintillation::sample_time ( G4double  tau1,
G4double  tau2 
)
private

Definition at line 711 of file G4Scintillation.cc.

712 {
713 // tau1: rise time and tau2: decay time
714 
715  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
716  while(1) {
717  // two random numbers
718  G4double ran1 = G4UniformRand();
719  G4double ran2 = G4UniformRand();
720  //
721  // exponential distribution as envelope function: very efficient
722  //
723  G4double d = (tau1+tau2)/tau2;
724  // make sure the envelope function is
725  // always larger than the bi-exponential
726  G4double t = -1.0*tau2*std::log(1-ran1);
727  G4double gg = d*single_exp(t,tau2);
728  if (ran2 <= bi_exp(t,tau1,tau2)/gg) return t;
729  }
730  return -1.0;
731 }
G4double bi_exp(G4double t, G4double tau1, G4double tau2)
Float_t d
G4double single_exp(G4double t, G4double tau2)
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool  state)

Definition at line 647 of file G4Scintillation.cc.

648 {
649  fFiniteRiseTime = state;
650 }
Here is the caller graph for this function:

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 664 of file G4Scintillation.cc.

665 {
666  if (fEmSaturation) {
667  G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
668  JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
670  }
671  fScintillationByParticleType = scintType;
672 }
G4EmSaturation * fEmSaturation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool fScintillationByParticleType
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetScintillationExcitationRatio()

void G4Scintillation::SetScintillationExcitationRatio ( const G4double  ratio)

Definition at line 657 of file G4Scintillation.cc.

657  {
658  fExcitationRatio = excitationratio;
659 }
G4double fExcitationRatio
Here is the caller graph for this function:

◆ SetScintillationYieldFactor()

void G4Scintillation::SetScintillationYieldFactor ( const G4double  yieldfactor)

Definition at line 652 of file G4Scintillation.cc.

653 {
654  fYieldFactor = yieldfactor;
655 }
Here is the caller graph for this function:

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool  state)

Definition at line 642 of file G4Scintillation.cc.

643 {
644  fTrackSecondariesFirst = state;
645 }
G4bool fTrackSecondariesFirst
Here is the caller graph for this function:

◆ single_exp()

G4double G4Scintillation::single_exp ( G4double  t,
G4double  tau2 
)
inlineprivate

Definition at line 327 of file G4Scintillation.hh.

328 {
329  return std::exp(-1.0*t/tau2)/tau2;
330 }
Here is the caller graph for this function:

Member Data Documentation

◆ fEmSaturation

G4EmSaturation* G4Scintillation::fEmSaturation
private

Definition at line 247 of file G4Scintillation.hh.

◆ fExcitationRatio

G4double G4Scintillation::fExcitationRatio
private

Definition at line 233 of file G4Scintillation.hh.

◆ fFastIntegralTable

G4PhysicsTable* G4Scintillation::fFastIntegralTable
protected

Definition at line 223 of file G4Scintillation.hh.

◆ fFiniteRiseTime

G4bool G4Scintillation::fFiniteRiseTime
private

Definition at line 229 of file G4Scintillation.hh.

◆ fScintillationByParticleType

G4bool G4Scintillation::fScintillationByParticleType
private

Definition at line 235 of file G4Scintillation.hh.

◆ fSlowIntegralTable

G4PhysicsTable* G4Scintillation::fSlowIntegralTable
protected

Definition at line 224 of file G4Scintillation.hh.

◆ fTrackSecondariesFirst

G4bool G4Scintillation::fTrackSecondariesFirst
private

Definition at line 228 of file G4Scintillation.hh.

◆ fYieldFactor

G4double G4Scintillation::fYieldFactor
private

Definition at line 231 of file G4Scintillation.hh.


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