Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // G4VProcess 27 // 28 // Class description: 29 // 30 // This class is the virtual class for physics process objects. 31 // It defines public methods which describe the behavior of 32 // a physics process. 33 34 // Authors: 35 // - 2 December 1995, G.Cosmo - First implementation, based on object model 36 // - 18 December 1996, H.Kurashige - New Physics scheme 37 // -------------------------------------------------------------------- 38 #ifndef G4VProcess_hh 39 #define G4VProcess_hh 1 40 41 #include <cmath> 42 43 #include "globals.hh" 44 #include "G4ios.hh" 45 #include "Randomize.hh" 46 47 #include "G4PhysicsTable.hh" 48 #include "G4VParticleChange.hh" 49 #include "G4ForceCondition.hh" 50 #include "G4GPILSelection.hh" 51 #include "G4ParticleChange.hh" 52 #include "G4ProcessType.hh" 53 54 class G4ParticleDefinition; 55 class G4DynamicParticle; 56 class G4Track; 57 class G4Step; 58 class G4ProcessTable; 59 60 class G4VProcess 61 { 62 63 public: 64 65 G4VProcess(const G4String& aName = "NoName", 66 G4ProcessType aType = fNotDefined); 67 // Constructor requires the process name and type 68 69 G4VProcess(const G4VProcess& right); 70 // Copy constructor copies the name but does not copy the 71 // physics table (null pointer is assigned instead) 72 73 virtual ~G4VProcess(); 74 // Destructor 75 76 G4VProcess& operator=(const G4VProcess&) = delete; 77 78 G4bool operator==(const G4VProcess& right) const; 79 G4bool operator!=(const G4VProcess& right) const; 80 // Equality operators 81 82 //////////////////////////// 83 // DoIt ///////////////// 84 //////////////////////////// 85 86 virtual G4VParticleChange* PostStepDoIt( 87 const G4Track& track, 88 const G4Step& stepData 89 ) = 0; 90 91 virtual G4VParticleChange* AlongStepDoIt( 92 const G4Track& track, 93 const G4Step& stepData 94 ) = 0; 95 virtual G4VParticleChange* AtRestDoIt( 96 const G4Track& track, 97 const G4Step& stepData 98 ) = 0; 99 // A virtual base class function that has to be overridden 100 // by any subclass. The DoIt() method actually performs the 101 // physics process and determines either momentum change 102 // of the production of secondaries etc. 103 // Arguments 104 // const G4Track& track: 105 // reference to the current G4Track information 106 // const G4Step& stepData: 107 // reference to the current G4Step information 108 109 ////////////////////////// 110 // GPIL /////////////// 111 ////////////////////////// 112 113 virtual G4double AlongStepGetPhysicalInteractionLength( 114 const G4Track& track, 115 G4double previousStepSize, 116 G4double currentMinimumStep, 117 G4double& proposedSafety, 118 G4GPILSelection* selection) = 0; 119 120 virtual G4double AtRestGetPhysicalInteractionLength( 121 const G4Track& track, 122 G4ForceCondition* condition ) = 0; 123 124 virtual G4double PostStepGetPhysicalInteractionLength( 125 const G4Track& track, 126 G4double previousStepSize, 127 G4ForceCondition* condition ) = 0; 128 // Returns the Step-size (actual length) which is allowed 129 // by "this" process. (for AtRestGetPhysicalInteractionLength, 130 // return value is Step-time) The NumberOfInteractionLengthLeft is 131 // recalculated by using previousStepSize and the Step-size is 132 // calucalted accoding to the resultant NumberOfInteractionLengthLeft. 133 // using NumberOfInteractionLengthLeft, which is recalculated at 134 // arguments 135 // const G4Track& track: 136 // reference to the current G4Track information 137 // G4double* previousStepSize: 138 // the Step-size (actual length) of the previous Step 139 // of this track. Negative calue indicates that 140 // NumberOfInteractionLengthLeft must be reset. 141 // the current physical interaction legth of this process 142 // G4ForceCondition* condition: 143 // the flag indicates DoIt of this process is forced 144 // to be called 145 // Forced: Corresponding DoIt is forced 146 // NotForced: Corresponding DoIt is called 147 // if the Step size of this Step is determined 148 // by this process 149 // !! AlongStepDoIt is always called !! 150 // G4double& currentMinimumStep: 151 // this value is used for transformation of 152 // true path length to geometrical path length 153 154 inline G4double GetCurrentInteractionLength() const; 155 // Returns currentInteractionLength 156 157 ////////// PIL factor //////// 158 // 159 inline void SetPILfactor(G4double value); 160 inline G4double GetPILfactor() const; 161 // Set/Get factor for PhysicsInteractionLength 162 // which is passed to G4SteppingManager for both AtRest and PostStep 163 164 // These three GPIL methods are used by Stepping Manager. 165 // They invoke virtual GPIL methods listed above. 166 // As for AtRest and PostStep the returned value is multipled by 167 // thePILfactor 168 // 169 inline G4double AlongStepGPIL( const G4Track& track, 170 G4double previousStepSize, 171 G4double currentMinimumStep, 172 G4double& proposedSafety, 173 G4GPILSelection* selection ); 174 175 inline G4double AtRestGPIL( const G4Track& track, 176 G4ForceCondition* condition ); 177 178 inline G4double PostStepGPIL( const G4Track& track, 179 G4double previousStepSize, 180 G4ForceCondition* condition ); 181 182 virtual G4bool IsApplicable(const G4ParticleDefinition&) { return true; } 183 // Returns true if this process object is applicable to 184 // the particle type. Process will not be registered to a 185 // particle if IsApplicable is false 186 187 virtual void BuildPhysicsTable(const G4ParticleDefinition&) {} 188 // Messaged by the Particle definition (via the Process manager) 189 // whenever cross-section tables have to be rebuilt (i.e. if new 190 // materials have been defined). 191 // It is overloaded by individual processes when they need physics 192 // tables 193 194 virtual void PreparePhysicsTable(const G4ParticleDefinition&) {} 195 // Messaged by the Particle definition (via the Process manager) 196 // whenever cross-section tables have to be prepared for rebuild 197 // (i.e. if new materials have been defined). 198 // It is overloaded by individual processes when they need physics 199 // tables 200 201 // Processes which Build physics tables independent of cuts 202 // (for example in their constructors) should preferably use private 203 // void BuildThePhysicsTable() and void PreparePhysicsTable(). 204 // *Not* another BuildPhysicsTable 205 206 virtual G4bool StorePhysicsTable(const G4ParticleDefinition* , 207 const G4String&, G4bool) { return true; } 208 // Store PhysicsTable in a file. 209 // Return false in case of failure at I/O 210 211 virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition* , 212 const G4String&, G4bool) { return false; } 213 // Retrieve Physics from a file. 214 // Return true if the Physics Table can be built by using file. 215 // Return false if the process has no functionality or in case 216 // of failure. File name should be defined by each process and the 217 // file should be placed under the directory specified by the argument 218 219 const G4String& GetPhysicsTableFileName(const G4ParticleDefinition* , 220 const G4String& directory, 221 const G4String& tableName, 222 G4bool ascii = false); 223 // This method is utility for Store/RetreivePhysicsTable 224 225 inline const G4String& GetProcessName() const; 226 // Returns the name of the process 227 228 inline G4ProcessType GetProcessType() const; 229 // Returns the process type 230 231 inline void SetProcessType(G4ProcessType); 232 // Sets the process type 233 234 inline G4int GetProcessSubType() const; 235 // Returns the process sub type 236 237 inline void SetProcessSubType(G4int); 238 // Sets the process sub type 239 240 static const G4String& GetProcessTypeName(G4ProcessType); 241 // Returns the process type name 242 243 virtual const G4VProcess* GetCreatorProcess() const; 244 // Returns the process to be used as CreatorProcess for secondaries 245 // coming from this process 246 247 virtual void StartTracking(G4Track*); 248 virtual void EndTracking(); 249 // Inform Start/End of tracking for each track to the physics process 250 251 virtual void SetProcessManager(const G4ProcessManager*); 252 // A process manager sets its own pointer when the process 253 // is registered in the process Manager 254 virtual const G4ProcessManager* GetProcessManager(); 255 // Get the process manager which the process belongs to 256 257 virtual void ResetNumberOfInteractionLengthLeft(); 258 // Reset (determine the value of) NumberOfInteractionLengthLeft 259 260 inline G4double GetNumberOfInteractionLengthLeft() const; 261 // Get NumberOfInteractionLengthLeft 262 263 inline G4double GetTotalNumberOfInteractionLengthTraversed() const; 264 // Get NumberOfInteractionLength after 265 // ResetNumberOfInteractionLengthLeft() is invoked 266 267 inline G4bool isAtRestDoItIsEnabled() const; 268 inline G4bool isAlongStepDoItIsEnabled() const; 269 inline G4bool isPostStepDoItIsEnabled() const; 270 // These methods indicate which DoIt is enabled. 271 // They are used by G4ProcessManager to check 272 // that ordering parameters are properly set 273 274 virtual void DumpInfo() const; 275 // Dump out process information 276 277 virtual void ProcessDescription(std::ostream& outfile) const; 278 // Write out to html file for automatic documentation 279 280 inline void SetVerboseLevel(G4int value); 281 inline G4int GetVerboseLevel() const; 282 // set/get control flag for output message 283 // 0: Silent 284 // 1: Warning message 285 // 2: More 286 287 virtual void SetMasterProcess(G4VProcess* masterP); 288 // Sets the master thread process instance 289 inline const G4VProcess* GetMasterProcess() const; 290 // Returns the master thread process instance. 291 // Can be used to initialise worker type processes 292 // instances from master one (e.g. to share a read-only table) 293 // if ( this != GetMasterProcess() ) { /*worker*/ } 294 // else { /* master or sequential */ } 295 296 virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& part); 297 // Messaged by the Particle definition (via the Process manager) 298 // in worker threads. See BuildWorkerPhyiscsTable() method. 299 // Can be used to share among threads physics tables. 300 // Use GetMasterProcess() to get pointer of master process from 301 // worker thread. 302 // By default this method makes a forward call to BuildPhysicsTable() 303 304 virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition&); 305 // Messaged by the Particle definition (via the Process manager) 306 // in worker threads. See PreparephysicsTable(). 307 // Can be used to share among threads physics tables. 308 // Use GetMasterProcess() to get pointer of master process from 309 // worker thread 310 // By default this method makes a forward call to PreparePhysicsTable() 311 312 protected: 313 314 inline void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize); 315 // Subtract NumberOfInteractionLengthLeft by the value corresponding 316 // to previousStepSize 317 318 inline void ClearNumberOfInteractionLengthLeft(); 319 // This method should be at the end of PostStepDoIt() and AtRestDoIt()! 320 321 protected: 322 323 const G4ProcessManager* aProcessManager = nullptr; 324 325 G4VParticleChange* pParticleChange = nullptr; 326 // The pointer to G4VParticleChange object 327 // which is modified and returned by address by the DoIt() method. 328 // This pointer should be set in each physics process 329 // after construction of derived class object 330 331 G4ParticleChange aParticleChange; 332 // This object is kept for compatibility with old scheme. 333 // May be removed in future 334 335 G4double theNumberOfInteractionLengthLeft = -1.0; 336 // The flight length left for the current tracking particle 337 // in unit of "Interaction length" 338 339 G4double currentInteractionLength = -1.0; 340 // The InteractionLength in the current material 341 342 G4double theInitialNumberOfInteractionLength = -1.0; 343 // The initial value when ResetNumberOfInteractionLengthLeft() is invoked 344 345 G4String theProcessName; 346 // The name of the process 347 348 G4String thePhysicsTableFileName; 349 350 G4ProcessType theProcessType = fNotDefined; 351 // The type of the process 352 353 G4int theProcessSubType = -1; 354 // The sub type of the process 355 356 G4double thePILfactor = 1.0; 357 // Factor for PhysicsInteractionLength 358 // which is passed to G4SteppingManager 359 360 G4int verboseLevel = 0; 361 // Controle flag for output message 362 363 G4bool enableAtRestDoIt = true; 364 G4bool enableAlongStepDoIt = true; 365 G4bool enablePostStepDoIt = true; 366 367 private: 368 369 G4VProcess(); 370 // Hidden default constructor 371 372 private: 373 374 G4VProcess* masterProcessShadow = nullptr; 375 // For multi-threaded: pointer to the instance of this process 376 // for the master thread 377 378 G4ProcessTable* fProcessTable = nullptr; 379 }; 380 381 // ----------------------------------------- 382 // inlined function members implementation 383 // ----------------------------------------- 384 385 inline 386 const G4String& G4VProcess::GetProcessName() const 387 { 388 return theProcessName; 389 } 390 391 inline 392 G4ProcessType G4VProcess::GetProcessType() const 393 { 394 return theProcessType; 395 } 396 397 inline 398 void G4VProcess::SetProcessType(G4ProcessType aType) 399 { 400 theProcessType = aType; 401 } 402 403 inline 404 G4int G4VProcess::GetProcessSubType() const 405 { 406 return theProcessSubType; 407 } 408 409 inline 410 void G4VProcess::SetProcessSubType(G4int value) 411 { 412 theProcessSubType = value; 413 } 414 415 inline 416 void G4VProcess::SetVerboseLevel(G4int value) 417 { 418 verboseLevel = value; 419 } 420 421 inline 422 G4int G4VProcess::GetVerboseLevel() const 423 { 424 return verboseLevel; 425 } 426 427 inline 428 void G4VProcess::ClearNumberOfInteractionLengthLeft() 429 { 430 theInitialNumberOfInteractionLength = -1.0; 431 theNumberOfInteractionLengthLeft = -1.0; 432 } 433 434 inline 435 G4double G4VProcess::GetNumberOfInteractionLengthLeft() const 436 { 437 return theNumberOfInteractionLengthLeft; 438 } 439 440 inline 441 G4double G4VProcess::GetTotalNumberOfInteractionLengthTraversed() const 442 { 443 return theInitialNumberOfInteractionLength - theNumberOfInteractionLengthLeft; 444 } 445 446 inline 447 G4double G4VProcess::GetCurrentInteractionLength() const 448 { 449 return currentInteractionLength; 450 } 451 452 inline 453 void G4VProcess::SetPILfactor(G4double value) 454 { 455 if (value>0.) { thePILfactor = value; } 456 } 457 458 inline 459 G4double G4VProcess::GetPILfactor() const 460 { 461 return thePILfactor; 462 } 463 464 inline 465 G4double G4VProcess::AlongStepGPIL( const G4Track& track, 466 G4double previousStepSize, 467 G4double currentMinimumStep, 468 G4double& proposedSafety, 469 G4GPILSelection* selection ) 470 { 471 return AlongStepGetPhysicalInteractionLength(track, previousStepSize, 472 currentMinimumStep, proposedSafety, selection); 473 } 474 475 inline 476 G4double G4VProcess::AtRestGPIL( const G4Track& track, 477 G4ForceCondition* condition ) 478 { 479 return thePILfactor * AtRestGetPhysicalInteractionLength(track, condition); 480 } 481 482 inline 483 G4double G4VProcess::PostStepGPIL( const G4Track& track, 484 G4double previousStepSize, 485 G4ForceCondition* condition ) 486 { 487 return thePILfactor * 488 PostStepGetPhysicalInteractionLength(track, previousStepSize, condition); 489 } 490 491 inline 492 void G4VProcess::SetProcessManager(const G4ProcessManager* procMan) 493 { 494 aProcessManager = procMan; 495 } 496 497 inline 498 const G4ProcessManager* G4VProcess::GetProcessManager() 499 { 500 return aProcessManager; 501 } 502 503 inline 504 G4bool G4VProcess::isAtRestDoItIsEnabled() const 505 { 506 return enableAtRestDoIt; 507 } 508 509 inline 510 G4bool G4VProcess::isAlongStepDoItIsEnabled() const 511 { 512 return enableAlongStepDoIt; 513 } 514 515 inline 516 G4bool G4VProcess::isPostStepDoItIsEnabled() const 517 { 518 return enablePostStepDoIt; 519 } 520 521 inline 522 const G4VProcess* G4VProcess::GetMasterProcess() const 523 { 524 return masterProcessShadow; 525 } 526 527 inline 528 void G4VProcess::SubtractNumberOfInteractionLengthLeft( G4double prevStepSize ) 529 { 530 if (currentInteractionLength>0.0) 531 { 532 theNumberOfInteractionLengthLeft -= prevStepSize/currentInteractionLength; 533 if(theNumberOfInteractionLengthLeft<0.) 534 { 535 theNumberOfInteractionLengthLeft=CLHEP::perMillion; 536 } 537 } 538 else 539 { 540 #ifdef G4VERBOSE 541 if (verboseLevel>0) 542 { 543 G4cerr << "G4VProcess::SubtractNumberOfInteractionLengthLeft()"; 544 G4cerr << " [" << theProcessName << "]" <<G4endl; 545 G4cerr << " currentInteractionLength = " 546 << currentInteractionLength << " [mm]"; 547 G4cerr << " previousStepSize = " << prevStepSize << " [mm]"; 548 G4cerr << G4endl; 549 } 550 #endif 551 G4String msg = "Negative currentInteractionLength for "; 552 msg += theProcessName; 553 G4Exception("G4VProcess::SubtractNumberOfInteractionLengthLeft()", 554 "ProcMan201", EventMustBeAborted, msg); 555 } 556 } 557 558 #endif 559