File indexing completed on 2026-04-10 07:50:30
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069 #ifndef Local_DsG4Scintillation_h
0070 #define Local_DsG4Scintillation_h 1
0071
0072
0073 #ifdef WITH_G4OPTICKS
0074 #include "plog/Severity.h"
0075 #endif
0076
0077 #include "globals.hh"
0078 #include "templates.hh"
0079 #include "Randomize.hh"
0080 #include "G4Poisson.hh"
0081 #include "G4ThreeVector.hh"
0082 #include "G4ParticleMomentum.hh"
0083 #include "G4Step.hh"
0084 #include "G4VRestDiscreteProcess.hh"
0085 #include "G4OpticalPhoton.hh"
0086 #include "G4DynamicParticle.hh"
0087 #include "G4Material.hh"
0088 #include "G4PhysicsTable.hh"
0089 #include "G4MaterialPropertiesTable.hh"
0090 #include "G4PhysicsOrderedFreeVector.hh"
0091 #include "G4UImessenger.hh"
0092
0093 #ifdef STANDALONE
0094 #else
0095 #include "DsPhysConsOptical.h"
0096 #endif
0097
0098 class G4UIcommand;
0099 class G4UIdirectory;
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109
0110
0111 class Local_DsG4Scintillation : public G4VRestDiscreteProcess, public G4UImessenger
0112 {
0113
0114 private:
0115 #ifdef WITH_G4OPTICKS
0116 static const plog::Severity LEVEL ;
0117 #endif
0118
0119 public:
0120 static const bool FLOAT ;
0121 static const int PIDX ;
0122 void ResetNumberOfInteractionLengthLeft();
0123
0124
0125
0126
0127
0128
0129
0130
0131 public:
0132
0133
0134
0135
0136
0137 Local_DsG4Scintillation(G4int opticksMode = 0, const G4String& processName = "Scintillation",
0138 G4ProcessType type = fElectromagnetic);
0139
0140
0141
0142 ~Local_DsG4Scintillation();
0143
0144
0145
0146
0147
0148 public:
0149
0150
0151
0152
0153
0154 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
0155
0156
0157
0158 G4double GetMeanFreePath(const G4Track& aTrack,
0159 G4double ,
0160 G4ForceCondition* );
0161
0162
0163
0164
0165 G4double GetMeanLifeTime(const G4Track& aTrack,
0166 G4ForceCondition* );
0167
0168
0169
0170
0171 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
0172 const G4Step& aStep);
0173 G4VParticleChange* AtRestDoIt (const G4Track& aTrack,
0174 const G4Step& aStep);
0175
0176
0177
0178 void SetTrackSecondariesFirst(const G4bool state);
0179
0180
0181
0182
0183 G4bool GetTrackSecondariesFirst() const;
0184
0185
0186
0187 void SetScintillationYieldFactor(const G4double yieldfactor);
0188
0189
0190
0191 G4double GetScintillationYieldFactor() const;
0192
0193
0194 void SetUseFastMu300nsTrick(const G4bool fastMu300nsTrick);
0195 G4bool GetUseFastMu300nsTrick() const;
0196
0197 #ifdef WITH_G4OPTICKS
0198 G4MaterialPropertyVector* getMaterialProperty(const char* name, G4int materialIndex) ;
0199 G4PhysicsOrderedFreeVector* getScintillationIntegral(G4int scnt, G4int materialIndex) const;
0200 G4double getSampledEnergy(G4int scnt, G4int materialIndex) const ;
0201 G4double getSampledWavelength(G4int scnt, G4int materialIndex) const ;
0202 #endif
0203
0204 void SetScintillationExcitationRatio(const G4double excitationratio);
0205
0206
0207
0208
0209
0210 G4double GetScintillationExcitationRatio() const;
0211
0212
0213 G4PhysicsTable* GetFastIntegralTable() const;
0214
0215
0216 G4PhysicsTable* GetSlowIntegralTable() const;
0217
0218
0219 G4PhysicsTable* GetReemissionIntegralTable() const;
0220
0221
0222 void DumpPhysicsTable() const;
0223
0224
0225
0226 G4double GetPhotonWeight() const { return fPhotonWeight; }
0227 void SetPhotonWeight(G4double weight) { fPhotonWeight = weight; }
0228 void SetDoReemission(bool tf = true) { doReemission = tf; }
0229 bool GetDoReemission() { return doReemission; }
0230 void SetDoBothProcess(bool tf = true) { doBothProcess = tf; }
0231 bool GetDoBothProcess() { return doBothProcess; }
0232
0233
0234 void SetDoReemissionOnly(bool tf = true) { doReemissionOnly = tf; }
0235 bool GetDoReemissionOnly() { return doReemissionOnly; }
0236
0237
0238 void SetDoQuenching(const G4bool enable) { fEnableQuenching = enable; }
0239 G4bool GetDoQuenching() const { return fEnableQuenching; }
0240
0241 G4bool GetApplyPreQE() const { return fApplyPreQE; }
0242 void SetApplyPreQE(G4bool a) { fApplyPreQE = a; }
0243
0244 G4double GetPreQE() const { return fPreQE; }
0245 void SetPreQE(G4double a) { fPreQE = a; }
0246
0247 void SetBirksConstant1(double c1) {birksConstant1 = c1;}
0248 double GetBirksConstant1() {return birksConstant1;}
0249 void SetBirksConstant2(double c2) {birksConstant2 = c2;}
0250 double GetBirksConstant2() {return birksConstant2;}
0251
0252 void SetSlowerTimeConstant(double st) {slowerTimeConstant = st;}
0253 double GetSlowerTimeConstant() {return slowerTimeConstant;}
0254
0255 void SetSlowerRatio(double sr) {slowerRatio = sr;}
0256 double GetSlowerRatio() {return slowerRatio;}
0257
0258
0259 void SetGammaSlowerTimeConstant(double st) { gammaSlowerTime = st;}
0260 double GetGammaSlowerTimeConstant() {return gammaSlowerTime;}
0261
0262 void SetGammaSlowerRatio(double sr) { gammaSlowerRatio = sr;}
0263 double GetGammaSlowerRatio() {return gammaSlowerRatio;}
0264
0265
0266 void SetNeutronSlowerTimeConstant(double st) { neutronSlowerTime = st;}
0267 double GetNeutronSlowerTimeConstant() {return neutronSlowerTime;}
0268
0269 void SetNeutronSlowerRatio(double sr) { neutronSlowerRatio = sr;}
0270 double GetNeutronSlowerRatio() {return neutronSlowerRatio;}
0271
0272
0273 void SetAlphaSlowerTimeConstant(double st) { alphaSlowerTime = st;}
0274 double GetAlphaSlowerTimeConstant() {return alphaSlowerTime;}
0275
0276 void SetAlphaSlowerRatio(double sr) { alphaSlowerRatio = sr;}
0277 double GetAlphaSlowerRatio() {return alphaSlowerRatio;}
0278
0279 void SetFlagDecayTimeFast(bool flag) { flagDecayTimeFast = flag; }
0280 bool GetFlagDecayTimeFast() { return flagDecayTimeFast; }
0281 void SetFlagDecayTimeSlow(bool flag) { flagDecayTimeSlow = flag; }
0282 bool GetFlagDecayTimeSlow() { return flagDecayTimeSlow; }
0283
0284
0285
0286
0287
0288 void SetNoOp(bool tf = true) { m_noop = tf; }
0289
0290 public:
0291
0292 G4PhysicsTable* getSlowIntegralTable();
0293 G4PhysicsTable* getFastIntegralTable();
0294 G4PhysicsTable* getReemissionIntegralTable();
0295
0296 private:
0297
0298 void BuildThePhysicsTable();
0299
0300
0301
0302
0303
0304
0305
0306 protected:
0307
0308 G4PhysicsTable* theSlowIntegralTable;
0309 G4PhysicsTable* theFastIntegralTable;
0310 G4PhysicsTable* theReemissionIntegralTable;
0311
0312
0313 G4bool doReemission;
0314
0315 G4bool doBothProcess;
0316
0317
0318
0319 G4bool doReemissionOnly;
0320
0321
0322 G4bool fEnableQuenching;
0323
0324
0325 double birksConstant1;
0326 double birksConstant2;
0327
0328 double slowerTimeConstant;
0329 double slowerRatio;
0330
0331 double gammaSlowerTime;
0332 double gammaSlowerRatio;
0333 double neutronSlowerTime;
0334 double neutronSlowerRatio;
0335 double alphaSlowerTime;
0336 double alphaSlowerRatio;
0337
0338
0339
0340
0341 bool flagDecayTimeFast;
0342 bool flagDecayTimeSlow;
0343
0344 private:
0345
0346 G4bool fTrackSecondariesFirst;
0347
0348 G4double YieldFactor;
0349 G4bool FastMu300nsTrick;
0350 G4double ExcitationRatio;
0351
0352
0353 G4double fPhotonWeight;
0354 G4bool fApplyPreQE;
0355 G4double fPreQE;
0356 bool m_noop;
0357 G4int m_opticksMode ;
0358 };
0359
0360
0361
0362
0363
0364
0365
0366
0367 inline
0368 G4bool Local_DsG4Scintillation::IsApplicable(const G4ParticleDefinition& aParticleType)
0369 {
0370 if (aParticleType.GetParticleName() == "opticalphoton"){
0371 return true;
0372 } else if (doReemissionOnly) {
0373
0374 return false;
0375 } else {
0376 return true;
0377 }
0378 }
0379
0380 inline
0381 void Local_DsG4Scintillation::SetTrackSecondariesFirst(const G4bool state)
0382 {
0383 fTrackSecondariesFirst = state;
0384 }
0385
0386 inline
0387 G4bool Local_DsG4Scintillation::GetTrackSecondariesFirst() const
0388 {
0389 return fTrackSecondariesFirst;
0390 }
0391
0392 inline
0393 void Local_DsG4Scintillation::SetScintillationYieldFactor(const G4double yieldfactor)
0394 {
0395 YieldFactor = yieldfactor;
0396 }
0397
0398
0399 inline
0400 G4double Local_DsG4Scintillation::GetScintillationYieldFactor() const
0401 {
0402 return YieldFactor;
0403 }
0404
0405 inline
0406 void Local_DsG4Scintillation::SetUseFastMu300nsTrick(const G4bool fastMu300nsTrick)
0407 {
0408 FastMu300nsTrick = fastMu300nsTrick;
0409 }
0410
0411 inline
0412 G4bool Local_DsG4Scintillation::GetUseFastMu300nsTrick() const
0413 {
0414 return FastMu300nsTrick;
0415 }
0416
0417
0418
0419
0420 inline
0421 void Local_DsG4Scintillation::SetScintillationExcitationRatio(const G4double excitationratio)
0422 {
0423 ExcitationRatio = excitationratio;
0424 }
0425
0426 inline
0427 G4double Local_DsG4Scintillation::GetScintillationExcitationRatio() const
0428 {
0429 return ExcitationRatio;
0430 }
0431
0432 inline
0433 G4PhysicsTable* Local_DsG4Scintillation::GetSlowIntegralTable() const
0434 {
0435 return theSlowIntegralTable;
0436 }
0437
0438 inline
0439 G4PhysicsTable* Local_DsG4Scintillation::GetFastIntegralTable() const
0440 {
0441 return theFastIntegralTable;
0442 }
0443
0444 inline
0445 G4PhysicsTable* Local_DsG4Scintillation::GetReemissionIntegralTable() const
0446 {
0447 return theReemissionIntegralTable;
0448 }
0449
0450 inline
0451 void Local_DsG4Scintillation::DumpPhysicsTable() const
0452 {
0453 if (theFastIntegralTable) {
0454 G4int PhysicsTableSize = theFastIntegralTable->entries();
0455 G4PhysicsOrderedFreeVector *v;
0456
0457 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
0458 {
0459 v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
0460 v->DumpValues();
0461 }
0462 }
0463
0464 if (theSlowIntegralTable) {
0465 G4int PhysicsTableSize = theSlowIntegralTable->entries();
0466 G4PhysicsOrderedFreeVector *v;
0467
0468 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
0469 {
0470 v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
0471 v->DumpValues();
0472 }
0473 }
0474
0475 if (theReemissionIntegralTable) {
0476 G4int PhysicsTableSize = theReemissionIntegralTable->entries();
0477 G4PhysicsOrderedFreeVector *v;
0478
0479 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
0480 {
0481 v = (G4PhysicsOrderedFreeVector*)(*theReemissionIntegralTable)[i];
0482 v->DumpValues();
0483 }
0484 }
0485 }
0486
0487 #endif