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

#include <G4OpticalParameters.hh>

Public Member Functions

 ~G4OpticalParameters ()
 
void SetDefaults ()
 
void StreamInfo (std::ostream &os) const
 
void Dump () const
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel () const
 
void SetProcessActivation (const G4String &, G4bool)
 
G4bool GetProcessActivation (const G4String &) const
 
void Configure (G4OpticalProcessIndex, G4bool)
 
G4bool GetConfiguration (G4OpticalProcessIndex)
 
void SetTrackSecondariesFirst (G4OpticalProcessIndex, G4bool)
 
G4bool GetTrackSecondariesFirst (G4OpticalProcessIndex)
 
void SetCerenkovMaxPhotonsPerStep (G4int)
 
G4int GetCerenkovMaxPhotonsPerStep () const
 
void SetCerenkovVerboseLevel (G4int)
 
G4int GetCerenkovVerboseLevel () const
 
void SetCerenkovMaxBetaChange (G4double)
 
G4double GetCerenkovMaxBetaChange () const
 
void SetCerenkovTrackSecondariesFirst (G4bool)
 
G4bool GetCerenkovTrackSecondariesFirst () const
 
void SetCerenkovStackPhotons (G4bool)
 
G4bool GetCerenkovStackPhotons () const
 
void SetScintYieldFactor (G4double)
 
G4double GetScintYieldFactor () const
 
void SetScintExcitationRatio (G4double)
 
G4double GetScintExcitationRatio () const
 
void SetScintByParticleType (G4bool)
 
G4bool GetScintByParticleType () const
 
void SetScintTrackInfo (G4bool)
 
G4bool GetScintTrackInfo () const
 
void SetScintTrackSecondariesFirst (G4bool)
 
G4bool GetScintTrackSecondariesFirst () const
 
void SetScintFiniteRiseTime (G4bool)
 
G4bool GetScintFiniteRiseTime () const
 
void SetScintStackPhotons (G4bool)
 
G4bool GetScintStackPhotons () const
 
void SetScintVerboseLevel (G4int)
 
G4int GetScintVerboseLevel () const
 
void SetScintEnhancedTimeConstants (G4bool)
 
G4bool GetScintEnhancedTimeConstants () const
 
void SetWLSTimeProfile (const G4String &)
 
G4String GetWLSTimeProfile () const
 
void SetWLSVerboseLevel (G4int)
 
G4int GetWLSVerboseLevel () const
 
void SetWLS2TimeProfile (const G4String &)
 
G4String GetWLS2TimeProfile () const
 
void SetWLS2VerboseLevel (G4int)
 
G4int GetWLS2VerboseLevel () const
 
void SetBoundaryVerboseLevel (G4int)
 
G4int GetBoundaryVerboseLevel () const
 
void SetBoundaryInvokeSD (G4bool)
 
G4bool GetBoundaryInvokeSD () const
 
void SetAbsorptionVerboseLevel (G4int)
 
G4int GetAbsorptionVerboseLevel () const
 
void SetRayleighVerboseLevel (G4int)
 
G4int GetRayleighVerboseLevel () const
 
void SetMieVerboseLevel (G4int)
 
G4int GetMieVerboseLevel () const
 

Static Public Member Functions

static G4OpticalParametersInstance ()
 

Friends

std::ostream & operator<< (std::ostream &os, const G4OpticalParameters &)
 

Detailed Description

Definition at line 96 of file G4OpticalParameters.hh.

Constructor & Destructor Documentation

◆ ~G4OpticalParameters()

G4OpticalParameters::~G4OpticalParameters ( )

Definition at line 78 of file G4OpticalParameters.cc.

79{
80 delete theMessenger;
81}

Member Function Documentation

◆ Configure()

void G4OpticalParameters::Configure ( G4OpticalProcessIndex  index,
G4bool  val 
)

Definition at line 191 of file G4OpticalParameters.cc.

192{
193 // DEPRECATED. Use SetProcessActivation instead.
194
195 // Configure the physics constructor to use/not use a selected process.
196 // This method can only be called in PreInit> phase (before execution of
197 // ConstructProcess). The process is not added to particle's process manager
198 // and so it cannot be re-activated later in Idle> phase with the command
199 // /process/activate.
200
201 if(IsLocked()) { return; }
202 if (index == kCerenkov) processActivation["Cerenkov"] = val;
203 else if (index == kScintillation) processActivation["Scintillation"] = val;
204 else if (index == kAbsorption) processActivation["OpAbsorption"] = val;
205 else if (index == kRayleigh) processActivation["OpRayleigh"] = val;
206 else if (index == kMieHG) processActivation["OpMieHG"] = val;
207 else if (index == kWLS) processActivation["OpWLS"] = val;
208 else if (index == kWLS2) processActivation["OpWLS2"] = val;
209 else {
211 ed << "Process index " << index << " out of bounds.";
212 G4Exception("G4OpticalParameters::Configure()", "Optical010", FatalException, ed);
213 }
215 ed2 << "Method Configure(G4OpticalProcessIndex, G4bool) is deprecated "
216 << "and will be removed in a future Geant4 version. Please use "
217 << "SetProcessActivation(G4String, G4bool) instead.";
218 PrintWarning(ed2);
219}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ kWLS
Wave Length Shifting process index.
@ kScintillation
Scintillation process index.
@ kWLS2
Second Wave Length Shifting process index.
@ kRayleigh
Rayleigh scattering process index.
@ kAbsorption
Absorption process index.
@ kCerenkov
Cerenkov process index.
@ kMieHG
Mie scattering process index.

◆ Dump()

void G4OpticalParameters::Dump ( ) const

Definition at line 580 of file G4OpticalParameters.cc.

581{
582#ifdef G4MULTITHREADED
583 G4MUTEXLOCK(&opticalParametersMutex);
584#endif
586#ifdef G4MULTITHREADED
587 G4MUTEXUNLOCK(&opticalParametersMutex);
588#endif
589}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
G4GLOB_DLL std::ostream G4cout
void StreamInfo(std::ostream &os) const

Referenced by G4OpticalPhysics::PrintStatistics(), and G4OpticalParametersMessenger::SetNewValue().

◆ GetAbsorptionVerboseLevel()

G4int G4OpticalParameters::GetAbsorptionVerboseLevel ( ) const

Definition at line 514 of file G4OpticalParameters.cc.

515{
516 return absorptionVerboseLevel;
517}

◆ GetBoundaryInvokeSD()

G4bool G4OpticalParameters::GetBoundaryInvokeSD ( ) const

Definition at line 503 of file G4OpticalParameters.cc.

504{
505 return boundaryInvokeSD;
506}

Referenced by G4OpBoundaryProcess::Initialise().

◆ GetBoundaryVerboseLevel()

G4int G4OpticalParameters::GetBoundaryVerboseLevel ( ) const

Definition at line 492 of file G4OpticalParameters.cc.

493{
494 return boundaryVerboseLevel;
495}

Referenced by G4OpBoundaryProcess::Initialise().

◆ GetCerenkovMaxBetaChange()

G4double G4OpticalParameters::GetCerenkovMaxBetaChange ( ) const

Definition at line 327 of file G4OpticalParameters.cc.

328{
329 return cerenkovMaxBetaChange;
330}

Referenced by G4Cerenkov::Initialise().

◆ GetCerenkovMaxPhotonsPerStep()

G4int G4OpticalParameters::GetCerenkovMaxPhotonsPerStep ( ) const

Definition at line 316 of file G4OpticalParameters.cc.

317{
318 return cerenkovMaxPhotons;
319}

Referenced by G4Cerenkov::Initialise(), and G4OpticalParametersMessenger::SetNewValue().

◆ GetCerenkovStackPhotons()

G4bool G4OpticalParameters::GetCerenkovStackPhotons ( ) const

Definition at line 294 of file G4OpticalParameters.cc.

295{
296 return cerenkovStackPhotons;
297}

Referenced by G4Cerenkov::Initialise().

◆ GetCerenkovTrackSecondariesFirst()

G4bool G4OpticalParameters::GetCerenkovTrackSecondariesFirst ( ) const

Definition at line 338 of file G4OpticalParameters.cc.

339{
340 return cerenkovTrackSecondariesFirst;
341}

Referenced by G4Cerenkov::Initialise().

◆ GetCerenkovVerboseLevel()

G4int G4OpticalParameters::GetCerenkovVerboseLevel ( ) const

Definition at line 305 of file G4OpticalParameters.cc.

306{
307 return cerenkovVerboseLevel;
308}

Referenced by G4Cerenkov::Initialise().

◆ GetConfiguration()

G4bool G4OpticalParameters::GetConfiguration ( G4OpticalProcessIndex  index)

Definition at line 221 of file G4OpticalParameters.cc.

222{
223 // DEPRECATED. Use GetProcessActivation instead.
224 if (index == kCerenkov) return processActivation["Cerenkov"];
225 else if (index == kScintillation) return processActivation["Scintillation"];
226 else if (index == kAbsorption) return processActivation["OpAbsorption"];
227 else if (index == kRayleigh) return processActivation["OpRayleigh"];
228 else if (index == kMieHG) return processActivation["OpMieHG"];
229 else if (index == kWLS) return processActivation["OpWLS"];
230 else if (index == kWLS2) return processActivation["OpWLS2"];
231 else {
233 ed << "Process index " << index << " out of bounds.";
234 G4Exception("G4OpticalParameters::GetConfiguration()", "Optical011", JustWarning, ed);
235 }
237 ed2 << "Method GetConfiguration(G4OpticalProcessIndex) is deprecated "
238 << "and will be removed in a future Geant4 version. Please use "
239 << "GetProcessActivation(G4String) instead.";
240 PrintWarning(ed2);
241 return true;
242}
@ JustWarning

◆ GetMieVerboseLevel()

G4int G4OpticalParameters::GetMieVerboseLevel ( ) const

Definition at line 536 of file G4OpticalParameters.cc.

537{
538 return mieVerboseLevel;
539}

◆ GetProcessActivation()

G4bool G4OpticalParameters::GetProcessActivation ( const G4String process) const

Definition at line 188 of file G4OpticalParameters.cc.

189{ return processActivation.find(process)->second; }

Referenced by StreamInfo().

◆ GetRayleighVerboseLevel()

G4int G4OpticalParameters::GetRayleighVerboseLevel ( ) const

Definition at line 525 of file G4OpticalParameters.cc.

526{
527 return rayleighVerboseLevel;
528}

◆ GetScintByParticleType()

G4bool G4OpticalParameters::GetScintByParticleType ( ) const

Definition at line 371 of file G4OpticalParameters.cc.

372{
373 return scintByParticleType;
374}

Referenced by G4Scintillation::Initialise().

◆ GetScintEnhancedTimeConstants()

G4bool G4OpticalParameters::GetScintEnhancedTimeConstants ( ) const

Definition at line 437 of file G4OpticalParameters.cc.

438{
439 return scintEnhancedTimeConstants;
440}

Referenced by G4Scintillation::Initialise().

◆ GetScintExcitationRatio()

G4double G4OpticalParameters::GetScintExcitationRatio ( ) const

Definition at line 360 of file G4OpticalParameters.cc.

361{
362 return scintExcitationRatio;
363}

Referenced by G4Scintillation::Initialise().

◆ GetScintFiniteRiseTime()

G4bool G4OpticalParameters::GetScintFiniteRiseTime ( ) const

Definition at line 404 of file G4OpticalParameters.cc.

405{
406 return scintFiniteRiseTime;
407}

Referenced by G4Scintillation::Initialise().

◆ GetScintStackPhotons()

G4bool G4OpticalParameters::GetScintStackPhotons ( ) const

Definition at line 415 of file G4OpticalParameters.cc.

416{
417 return scintStackPhotons;
418}

Referenced by G4Scintillation::Initialise().

◆ GetScintTrackInfo()

G4bool G4OpticalParameters::GetScintTrackInfo ( ) const

Definition at line 382 of file G4OpticalParameters.cc.

383{
384 return scintTrackInfo;
385}

Referenced by G4Scintillation::Initialise().

◆ GetScintTrackSecondariesFirst()

G4bool G4OpticalParameters::GetScintTrackSecondariesFirst ( ) const

Definition at line 393 of file G4OpticalParameters.cc.

394{
395 return scintTrackSecondariesFirst;
396}

Referenced by G4Scintillation::Initialise().

◆ GetScintVerboseLevel()

G4int G4OpticalParameters::GetScintVerboseLevel ( ) const

Definition at line 426 of file G4OpticalParameters.cc.

427{
428 return scintVerboseLevel;
429}

Referenced by G4Scintillation::Initialise().

◆ GetScintYieldFactor()

G4double G4OpticalParameters::GetScintYieldFactor ( ) const

Definition at line 349 of file G4OpticalParameters.cc.

350{
351 return scintYieldFactor;
352}

Referenced by G4Scintillation::Initialise().

◆ GetTrackSecondariesFirst()

G4bool G4OpticalParameters::GetTrackSecondariesFirst ( G4OpticalProcessIndex  index)

Definition at line 267 of file G4OpticalParameters.cc.

270{
271 if (index == kCerenkov) return cerenkovTrackSecondariesFirst;
272 else if (index == kScintillation) return scintTrackSecondariesFirst;
273 else {
275 ed << "Process index " << index << " out of bounds.";
276 G4Exception("G4OpticalParameters::GetTrackSecondariesFirst()",
277 "Optical012", JustWarning, ed);
278 }
280 ed2 << "Method GetTrackSecondariesFirst(G4OpticalProcessIndex) is "
281 << "deprecated and will be removed in a future Geant4 version. Please use "
282 << "GetCerenkovTrackSecondariesFirst() and "
283 << "GetScintTrackSecondariesFirst() instead.";
284 PrintWarning(ed2);
285 return true;
286}

◆ GetVerboseLevel()

G4int G4OpticalParameters::GetVerboseLevel ( ) const

Definition at line 160 of file G4OpticalParameters.cc.

161{
162 return verboseLevel;
163}

◆ GetWLS2TimeProfile()

G4String G4OpticalParameters::GetWLS2TimeProfile ( ) const

Definition at line 470 of file G4OpticalParameters.cc.

471{
472 return wls2TimeProfileName;
473}

Referenced by G4OpWLS2::Initialise().

◆ GetWLS2VerboseLevel()

G4int G4OpticalParameters::GetWLS2VerboseLevel ( ) const

Definition at line 481 of file G4OpticalParameters.cc.

482{
483 return wls2VerboseLevel;
484}

Referenced by G4OpWLS2::Initialise().

◆ GetWLSTimeProfile()

G4String G4OpticalParameters::GetWLSTimeProfile ( ) const

Definition at line 448 of file G4OpticalParameters.cc.

449{
450 return wlsTimeProfileName;
451}

Referenced by G4OpWLS::Initialise().

◆ GetWLSVerboseLevel()

G4int G4OpticalParameters::GetWLSVerboseLevel ( ) const

Definition at line 459 of file G4OpticalParameters.cc.

460{
461 return wlsVerboseLevel;
462}

Referenced by G4OpWLS::Initialise().

◆ Instance()

G4OpticalParameters * G4OpticalParameters::Instance ( )
static

Definition at line 59 of file G4OpticalParameters.cc.

60{
61 if(nullptr == theInstance) {
62#ifdef G4MULTITHREADED
63 G4MUTEXLOCK(&opticalParametersMutex);
64 if(nullptr == theInstance) {
65#endif
66 static G4OpticalParameters manager;
67 theInstance = &manager;
68#ifdef G4MULTITHREADED
69 }
70 G4MUTEXUNLOCK(&opticalParametersMutex);
71#endif
72 }
73 return theInstance;
74}

Referenced by G4OpticalPhysics::Configure(), G4OpticalPhysics::ConstructProcess(), G4OpticalPhysics::G4OpticalPhysics(), G4Cerenkov::Initialise(), G4Scintillation::Initialise(), G4OpAbsorption::Initialise(), G4OpBoundaryProcess::Initialise(), G4OpMieHG::Initialise(), G4OpRayleigh::Initialise(), G4OpWLS::Initialise(), G4OpWLS2::Initialise(), G4OpticalPhysics::PrintStatistics(), G4OpticalPhysics::SetAbsorptionVerbosity(), G4OpticalPhysics::SetBoundaryVerbosity(), G4OpticalPhysics::SetCerenkovStackPhotons(), G4OpticalPhysics::SetCerenkovTrackSecondariesFirst(), G4OpticalPhysics::SetCerenkovVerbosity(), G4OpticalPhysics::SetFiniteRiseTime(), G4OpticalPhysics::SetInvokeSD(), G4OpticalPhysics::SetMaxBetaChangePerStep(), G4OpticalPhysics::SetMaxNumPhotonsPerStep(), G4OpticalPhysics::SetMieVerbosity(), G4OpticalPhysics::SetRayleighVerbosity(), G4OpticalPhysics::SetScintillationByParticleType(), G4OpticalPhysics::SetScintillationEnhancedTimeConstants(), G4OpticalPhysics::SetScintillationExcitationRatio(), G4OpticalPhysics::SetScintillationStackPhotons(), G4OpticalPhysics::SetScintillationTrackInfo(), G4OpticalPhysics::SetScintillationTrackSecondariesFirst(), G4OpticalPhysics::SetScintillationVerbosity(), G4OpticalPhysics::SetScintillationYieldFactor(), G4OpticalPhysics::SetTrackSecondariesFirst(), G4OpticalPhysics::SetWLSTimeProfile(), and G4OpticalPhysics::SetWLSVerbosity().

◆ SetAbsorptionVerboseLevel()

void G4OpticalParameters::SetAbsorptionVerboseLevel ( G4int  val)

Definition at line 508 of file G4OpticalParameters.cc.

509{
510 if(IsLocked()) {return;}
511 absorptionVerboseLevel = val;
512}

Referenced by G4OpticalPhysics::SetAbsorptionVerbosity(), G4OpticalParametersMessenger::SetNewValue(), and SetVerboseLevel().

◆ SetBoundaryInvokeSD()

void G4OpticalParameters::SetBoundaryInvokeSD ( G4bool  val)

Definition at line 497 of file G4OpticalParameters.cc.

498{
499 if(IsLocked()) {return;}
500 boundaryInvokeSD = val;
501}

Referenced by G4OpticalPhysics::SetInvokeSD(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetBoundaryVerboseLevel()

void G4OpticalParameters::SetBoundaryVerboseLevel ( G4int  val)

Definition at line 486 of file G4OpticalParameters.cc.

487{
488 if(IsLocked()) {return;}
489 boundaryVerboseLevel = val;
490}

Referenced by G4OpticalPhysics::SetBoundaryVerbosity(), G4OpticalParametersMessenger::SetNewValue(), and SetVerboseLevel().

◆ SetCerenkovMaxBetaChange()

void G4OpticalParameters::SetCerenkovMaxBetaChange ( G4double  val)

Definition at line 321 of file G4OpticalParameters.cc.

322{
323 if(IsLocked()) { return; }
324 cerenkovMaxBetaChange = val;
325}

Referenced by G4OpticalPhysics::SetMaxBetaChangePerStep(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetCerenkovMaxPhotonsPerStep()

void G4OpticalParameters::SetCerenkovMaxPhotonsPerStep ( G4int  val)

Definition at line 310 of file G4OpticalParameters.cc.

311{
312 if(IsLocked()) { return; }
313 cerenkovMaxPhotons = val;
314}

Referenced by G4OpticalPhysics::SetMaxNumPhotonsPerStep(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetCerenkovStackPhotons()

void G4OpticalParameters::SetCerenkovStackPhotons ( G4bool  val)

Definition at line 288 of file G4OpticalParameters.cc.

289{
290 if(IsLocked()) { return; }
291 cerenkovStackPhotons = val;
292}

Referenced by G4OpticalPhysics::SetCerenkovStackPhotons(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetCerenkovTrackSecondariesFirst()

void G4OpticalParameters::SetCerenkovTrackSecondariesFirst ( G4bool  val)

Definition at line 332 of file G4OpticalParameters.cc.

333{
334 if(IsLocked()) { return; }
335 cerenkovTrackSecondariesFirst = val;
336}

Referenced by G4OpticalPhysics::SetCerenkovTrackSecondariesFirst(), G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetTrackSecondariesFirst().

◆ SetCerenkovVerboseLevel()

void G4OpticalParameters::SetCerenkovVerboseLevel ( G4int  val)

Definition at line 299 of file G4OpticalParameters.cc.

300{
301 if(IsLocked()) { return; }
302 cerenkovVerboseLevel = val;
303}

Referenced by G4OpticalPhysics::SetCerenkovVerbosity(), G4OpticalParametersMessenger::SetNewValue(), and SetVerboseLevel().

◆ SetDefaults()

void G4OpticalParameters::SetDefaults ( )

Definition at line 93 of file G4OpticalParameters.cc.

94{
95 if(!IsLocked()) {
96 Initialise();
97 }
98}

◆ SetMieVerboseLevel()

void G4OpticalParameters::SetMieVerboseLevel ( G4int  val)

Definition at line 530 of file G4OpticalParameters.cc.

531{
532 if(IsLocked()) {return;}
533 mieVerboseLevel = val;
534}

Referenced by G4OpticalPhysics::SetMieVerbosity(), G4OpticalParametersMessenger::SetNewValue(), and SetVerboseLevel().

◆ SetProcessActivation()

void G4OpticalParameters::SetProcessActivation ( const G4String process,
G4bool  val 
)

Definition at line 165 of file G4OpticalParameters.cc.

166{
167 // Configure the physics constructor to use/not use a selected process.
168 // This method can only be called in PreInit> phase (before execution of
169 // ConstructProcess). The process is not added to particle's process manager
170 // and so it cannot be re-activated later in Idle> phase with the command
171 // /process/activate.
172
173 if(IsLocked()) { return; }
174 if (processActivation[process] == val) return;
175
176 // processActivation keys defined at initialisation
177 if (processActivation.find(process) != processActivation.end()) {
178 processActivation[process] = val;
179 }
180 else {
182 ed << "Process name " << process << " out of bounds.";
183 G4Exception("G4OpticalParameters::SetProcessActivation()", "Optical013",
184 FatalException, ed);
185 }
186}

Referenced by G4OpticalPhysics::Configure(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetRayleighVerboseLevel()

void G4OpticalParameters::SetRayleighVerboseLevel ( G4int  val)

Definition at line 519 of file G4OpticalParameters.cc.

520{
521 if(IsLocked()) {return;}
522 rayleighVerboseLevel = val;
523}

Referenced by G4OpticalParametersMessenger::SetNewValue(), G4OpticalPhysics::SetRayleighVerbosity(), and SetVerboseLevel().

◆ SetScintByParticleType()

void G4OpticalParameters::SetScintByParticleType ( G4bool  val)

Definition at line 365 of file G4OpticalParameters.cc.

366{
367 if(IsLocked()) { return; }
368 scintByParticleType = val;
369}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationByParticleType().

◆ SetScintEnhancedTimeConstants()

void G4OpticalParameters::SetScintEnhancedTimeConstants ( G4bool  val)

Definition at line 431 of file G4OpticalParameters.cc.

432{
433 if(IsLocked()) { return; }
434 scintEnhancedTimeConstants = val;
435}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationEnhancedTimeConstants().

◆ SetScintExcitationRatio()

void G4OpticalParameters::SetScintExcitationRatio ( G4double  val)

Definition at line 354 of file G4OpticalParameters.cc.

355{
356 if(IsLocked()) { return; }
357 scintExcitationRatio = val;
358}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationExcitationRatio().

◆ SetScintFiniteRiseTime()

void G4OpticalParameters::SetScintFiniteRiseTime ( G4bool  val)

Definition at line 398 of file G4OpticalParameters.cc.

399{
400 if(IsLocked()) { return; }
401 scintFiniteRiseTime = val;
402}

Referenced by G4OpticalPhysics::SetFiniteRiseTime(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetScintStackPhotons()

void G4OpticalParameters::SetScintStackPhotons ( G4bool  val)

Definition at line 409 of file G4OpticalParameters.cc.

410{
411 if(IsLocked()) { return; }
412 scintStackPhotons = val;
413}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationStackPhotons().

◆ SetScintTrackInfo()

void G4OpticalParameters::SetScintTrackInfo ( G4bool  val)

Definition at line 376 of file G4OpticalParameters.cc.

377{
378 if(IsLocked()) { return; }
379 scintTrackInfo = val;
380}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationTrackInfo().

◆ SetScintTrackSecondariesFirst()

void G4OpticalParameters::SetScintTrackSecondariesFirst ( G4bool  val)

Definition at line 387 of file G4OpticalParameters.cc.

388{
389 if(IsLocked()) { return; }
390 scintTrackSecondariesFirst = val;
391}

Referenced by G4OpticalParametersMessenger::SetNewValue(), G4OpticalPhysics::SetScintillationTrackSecondariesFirst(), and G4OpticalPhysics::SetTrackSecondariesFirst().

◆ SetScintVerboseLevel()

void G4OpticalParameters::SetScintVerboseLevel ( G4int  val)

Definition at line 420 of file G4OpticalParameters.cc.

421{
422 if(IsLocked()) { return; }
423 scintVerboseLevel = val;
424}

Referenced by G4OpticalParametersMessenger::SetNewValue(), G4OpticalPhysics::SetScintillationVerbosity(), and SetVerboseLevel().

◆ SetScintYieldFactor()

void G4OpticalParameters::SetScintYieldFactor ( G4double  val)

Definition at line 343 of file G4OpticalParameters.cc.

344{
345 if(IsLocked()) { return; }
346 scintYieldFactor = val;
347}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetScintillationYieldFactor().

◆ SetTrackSecondariesFirst()

void G4OpticalParameters::SetTrackSecondariesFirst ( G4OpticalProcessIndex  index,
G4bool  val 
)

Definition at line 244 of file G4OpticalParameters.cc.

246{
247 // DEPRECATED. Use SetCerenkovTrackSecondariesFirst and
248 // SetScintTrackSecondariesFirst instead.
249 if(IsLocked()) { return; }
250 if (index == kCerenkov) cerenkovTrackSecondariesFirst = val;
251 else if (index == kScintillation) scintTrackSecondariesFirst = val;
252 else {
254 ed << "Process index " << index << " out of bounds.";
255 G4Exception("G4OpticalParameters::SetTrackSecondariesFirst()",
256 "Optical013", FatalException, ed);
257 }
259 ed2 << "Method SetTrackSecondariesFirst(G4OpticalProcessIndex, G4bool) is "
260 << "deprecated and will be removed in a future Geant4 version. Please use "
261 << "SetCerenkovTrackSecondariesFirst(G4bool) and "
262 << "SetScintTrackSecondariesFirst(G4bool) instead.";
263 PrintWarning(ed2);
264
265}

◆ SetVerboseLevel()

void G4OpticalParameters::SetVerboseLevel ( G4int  val)

Definition at line 146 of file G4OpticalParameters.cc.

147{
148 if(IsLocked()) { return; }
149 verboseLevel = val;
150 SetCerenkovVerboseLevel(verboseLevel);
151 SetScintVerboseLevel(verboseLevel);
152 SetRayleighVerboseLevel(verboseLevel);
153 SetAbsorptionVerboseLevel(verboseLevel);
154 SetMieVerboseLevel(verboseLevel);
155 SetBoundaryVerboseLevel(verboseLevel);
156 SetWLSVerboseLevel(verboseLevel);
157 SetWLS2VerboseLevel(verboseLevel);
158}
void SetRayleighVerboseLevel(G4int)
void SetBoundaryVerboseLevel(G4int)
void SetAbsorptionVerboseLevel(G4int)
void SetCerenkovVerboseLevel(G4int)

Referenced by G4OpticalPhysics::G4OpticalPhysics(), and G4OpticalParametersMessenger::SetNewValue().

◆ SetWLS2TimeProfile()

void G4OpticalParameters::SetWLS2TimeProfile ( const G4String val)

Definition at line 464 of file G4OpticalParameters.cc.

465{
466 if(IsLocked()) { return; }
467 wls2TimeProfileName = val;
468}

Referenced by G4OpticalParametersMessenger::SetNewValue().

◆ SetWLS2VerboseLevel()

void G4OpticalParameters::SetWLS2VerboseLevel ( G4int  val)

Definition at line 475 of file G4OpticalParameters.cc.

476{
477 if(IsLocked()) {return; }
478 wls2VerboseLevel = val;
479}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and SetVerboseLevel().

◆ SetWLSTimeProfile()

void G4OpticalParameters::SetWLSTimeProfile ( const G4String val)

Definition at line 442 of file G4OpticalParameters.cc.

443{
444 if(IsLocked()) { return; }
445 wlsTimeProfileName = val;
446}

Referenced by G4OpticalParametersMessenger::SetNewValue(), and G4OpticalPhysics::SetWLSTimeProfile().

◆ SetWLSVerboseLevel()

void G4OpticalParameters::SetWLSVerboseLevel ( G4int  val)

Definition at line 453 of file G4OpticalParameters.cc.

454{
455 if(IsLocked()) {return; }
456 wlsVerboseLevel = val;
457}

Referenced by G4OpticalParametersMessenger::SetNewValue(), SetVerboseLevel(), and G4OpticalPhysics::SetWLSVerbosity().

◆ StreamInfo()

void G4OpticalParameters::StreamInfo ( std::ostream &  os) const

Definition at line 546 of file G4OpticalParameters.cc.

547{
548 G4int prec = os.precision(5);
549 os << "=======================================================================" << "\n";
550 os << "====== Optical Physics Parameters ========" << "\n";
551 os << "=======================================================================" << "\n";
552
553 os << " Cerenkov process active: " << GetProcessActivation("Cerenkov") << "\n";
554 os << " Cerenkov maximum photons per step: " << cerenkovMaxPhotons << "\n";
555 os << " Cerenkov maximum beta change per step: " << cerenkovMaxBetaChange << " %\n";
556 os << " Cerenkov stack photons: " << cerenkovStackPhotons << "\n";
557 os << " Cerenkov track secondaries first: " << cerenkovTrackSecondariesFirst << "\n";
558 os << " Scintillation process active: " << GetProcessActivation("Scintillation") << "\n";
559 os << " Scintillation yield factor: " << scintYieldFactor << "\n";
560 os << " Scintillation excitation ratio: " << scintExcitationRatio << "\n";
561 os << " Scintillation finite rise time: " << scintFiniteRiseTime << "\n";
562 os << " Scintillation by particle type: " << scintByParticleType << "\n";
563 os << " Scintillation record track info: " << scintTrackInfo << "\n";
564 os << " Scintillation stack photons: " << scintStackPhotons << "\n";
565 os << " Scintillation use enhanced time constants: " << scintEnhancedTimeConstants << "\n";
566 os << " Scintillation track secondaries first: " << scintTrackSecondariesFirst << "\n";
567 os << " WLS process active: " << GetProcessActivation("OpWLS") << "\n";
568 os << " WLS time profile name: " << wlsTimeProfileName << "\n";
569 os << " WLS2 process active: " << GetProcessActivation("OpWLS2") << "\n";
570 os << " WLS2 time profile name: " << wls2TimeProfileName << "\n";
571 os << " Boundary process active: " << GetProcessActivation("OpBoundary") << "\n";
572 os << " Boundary invoke sensitive detector: " << boundaryInvokeSD << "\n";
573 os << " Rayleigh process active: " << GetProcessActivation("OpRayleigh") << "\n";
574 os << " MieHG process active: " << GetProcessActivation("OpMieHG") << "\n";
575 os << " Absorption process active: " << GetProcessActivation("OpAbsorption") << "\n";
576 os << "=======================================================================" << "\n";
577 os.precision(prec);
578}
int G4int
Definition: G4Types.hh:85
G4bool GetProcessActivation(const G4String &) const

Referenced by Dump().

Friends And Related Function Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  os,
const G4OpticalParameters par 
)
friend

Definition at line 591 of file G4OpticalParameters.cc.

592{
593 par.StreamInfo(os);
594 return os;
595}

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