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 // 27 // 28 //--------------------------------------------------------------- 29 // 30 // G4FastSimulationManager.cc 31 // 32 // Description: 33 // Manages the Fast Simulation models attached to a envelope. 34 // 35 // History: 36 // Oct 97: Verderi && MoraDeFreitas - First Implementation. 37 // ... 38 // May 07: Move to parallel world scheme 39 // 40 //--------------------------------------------------------------- 41 42 #include "G4FastSimulationManager.hh" 43 44 #include "G4GlobalFastSimulationManager.hh" 45 #include "G4PVPlacement.hh" 46 #include "G4TransportationManager.hh" 47 48 // -------------------------------------------------- 49 // Constructor with envelope and IsUnique flag : 50 // -------------------------------------------------- 51 // 52 G4FastSimulationManager::G4FastSimulationManager(G4Envelope* anEnvelope, G4bool IsUnique) 53 : fFastTrack(anEnvelope, IsUnique) 54 { 55 // Communicates to the region that it becomes a 56 // envelope and with this fast simulation manager. 57 anEnvelope->SetFastSimulationManager(this); 58 59 // Add itself to the GlobalFastSimulationManager 60 G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->AddFastSimulationManager(this); 61 } 62 63 // ----------- 64 // Destructor: 65 // ----------- 66 G4FastSimulationManager::~G4FastSimulationManager() 67 { 68 // 69 // Check out the Envelope about this pointer. If in use, 70 // resets the Logical Volume IsEnvelope flag to avoid clash. 71 // 72 if (fFastTrack.GetEnvelope()->GetFastSimulationManager() == this) 73 fFastTrack.GetEnvelope()->ClearFastSimulationManager(); 74 // Remove itself from the GlobalFastSimulationManager 75 G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->RemoveFastSimulationManager( 76 this); 77 } 78 79 // --------------------------------------- 80 // Methods to activate/inactivate models 81 //---------------------------------------- 82 83 G4bool G4FastSimulationManager::ActivateFastSimulationModel(const G4String& aName) 84 { 85 G4int iModel; 86 87 // If the model is already active, do nothing. 88 for (iModel = 0; iModel < (G4int)ModelList.size(); ++iModel) 89 if (ModelList[iModel]->GetName() == aName) return true; 90 91 // Look for in the fInactivatedModels list, if found push_back it back to 92 // the ModelList 93 for (iModel = 0; iModel < (G4int)fInactivatedModels.size(); ++iModel) 94 if (fInactivatedModels[iModel]->GetName() == aName) { 95 ModelList.push_back(fInactivatedModels.removeAt(iModel)); 96 // forces the fApplicableModelList to be rebuild 97 fLastCrossedParticle = nullptr; 98 return true; 99 } 100 return false; 101 } 102 103 G4bool G4FastSimulationManager::InActivateFastSimulationModel(const G4String& aName) 104 { 105 // Look for in the ModelList, if found remove from it and keep the pointer 106 // on the fInactivatedModels list. 107 for (G4int iModel = 0; iModel < (G4int)ModelList.size(); ++iModel) 108 if (ModelList[iModel]->GetName() == aName) { 109 fInactivatedModels.push_back(ModelList.removeAt(iModel)); 110 // forces the fApplicableModelList to be rebuild 111 fLastCrossedParticle = nullptr; 112 return true; 113 } 114 return false; 115 } 116 117 G4VFastSimulationModel* 118 G4FastSimulationManager::GetFastSimulationModel(const G4String& modelName, 119 const G4VFastSimulationModel* previousFound, 120 G4bool& foundPrevious) const 121 { 122 G4VFastSimulationModel* model = nullptr; 123 for (auto iModel : ModelList) { 124 if (iModel->GetName() == modelName) { 125 if (previousFound == nullptr) { 126 model = iModel; 127 break; 128 } 129 if (iModel == previousFound) { 130 foundPrevious = true; 131 continue; 132 } 133 if (foundPrevious) { 134 model = iModel; 135 break; 136 } 137 } 138 } 139 return model; 140 } 141 142 void G4FastSimulationManager::FlushModels() 143 { 144 for (auto& iModel : ModelList) { 145 iModel->Flush(); 146 } 147 } 148 149 //------------------------------------------------------------------ 150 // Interface trigger method for the G4ParameterisationManagerProcess 151 //------------------------------------------------------------------ 152 // G4bool GetFastSimulationManagerTrigger(const G4Track &); 153 // 154 // This method is used to interface the G4FastSimulationManagerProcess 155 // with the user Fast Simulation Models. It's called when the particle 156 // is inside the envelope. 157 // 158 // It : 159 // 160 // 1) initialises the private members (fFastTrack and so 161 // on); 162 // 2) loops on the IsApplicable() methods to find out the 163 // ones should be applied. 164 // 2) for these, loops on the ModelTrigger() methods to find out 165 // perhaps one that must be applied just now. 166 // 167 // If the a Fast Simulation Model is triggered then it returns 168 // true, false otherwise. 169 // 170 //----------------------------------------------------------- 171 G4bool 172 G4FastSimulationManager::PostStepGetFastSimulationManagerTrigger(const G4Track& track, 173 const G4Navigator* theNavigator) 174 { 175 // If particle type changed re-build the fApplicableModelList. 176 if (fLastCrossedParticle != track.GetDefinition()) { 177 fLastCrossedParticle = track.GetDefinition(); 178 fApplicableModelList.clear(); 179 // If Model List is empty, do nothing ! 180 if (ModelList.empty()) return false; 181 182 for (auto iModel : ModelList) { 183 if (iModel->IsApplicable(*(track.GetDefinition()))) { 184 fApplicableModelList.push_back(iModel); 185 } 186 } 187 } 188 189 // If Applicable Model List is empty, do nothing ! 190 if (fApplicableModelList.empty()) return false; 191 192 // -- Register current track 193 fFastTrack.SetCurrentTrack(track, theNavigator); 194 195 // tests if particle are on the boundary and leaving, 196 // in this case do nothing ! 197 if (fFastTrack.OnTheBoundaryButExiting()) return false; 198 199 // Loops on the ModelTrigger() methods 200 for (auto iModel : fApplicableModelList) 201 //--------------------------------------------------- 202 // Asks the ModelTrigger method if it must be trigged now. 203 //--------------------------------------------------- 204 205 if (iModel->ModelTrigger(fFastTrack)) { 206 //-------------------------------------------------- 207 // The model will be applied. Initializes the G4FastStep 208 // with the current state of the G4Track and 209 // same usefull parameters. 210 // In particular it does SetLocalEnergyDeposit(0.0). 211 //-------------------------------------------------- 212 fFastStep.Initialize(fFastTrack); 213 214 // Keeps the FastSimulationModel pointer to call the 215 // DoIt() method. 216 fTriggedFastSimulationModel = iModel; 217 return true; 218 } 219 220 //-------------------------------------------- 221 // Nobody asks to gain control, returns false 222 //-------------------------------------------- 223 return false; 224 } 225 226 G4VParticleChange* G4FastSimulationManager::InvokePostStepDoIt() 227 { 228 fTriggedFastSimulationModel->DoIt(fFastTrack, fFastStep); 229 return &fFastStep; 230 } 231 232 // ------------------------------------------------------------- 233 // -- Mostly the same as above, in the case of AtRest particles: 234 // ------------------------------------------------------------- 235 G4bool 236 G4FastSimulationManager::AtRestGetFastSimulationManagerTrigger(const G4Track& track, 237 const G4Navigator* theNavigator) 238 { 239 // If particle type changed re-build the fApplicableModelList. 240 if (fLastCrossedParticle != track.GetDefinition()) { 241 fLastCrossedParticle = track.GetDefinition(); 242 fApplicableModelList.clear(); 243 // If Model List is empty, do nothing ! 244 if (ModelList.empty()) return false; 245 246 for (auto iModel : ModelList) { 247 if (iModel->IsApplicable(*(track.GetDefinition()))) { 248 fApplicableModelList.push_back(iModel); 249 } 250 } 251 } 252 253 // If Applicable Model List is empty, do nothing ! 254 if (fApplicableModelList.empty()) return false; 255 256 // -- Register current track 257 fFastTrack.SetCurrentTrack(track, theNavigator); 258 259 // -- (note: compared to the PostStepGetFastSimulationManagerTrigger, 260 // -- the test to see if the particle is on the boundary but leaving 261 // -- is irrelevant here) 262 263 // Loops on the models to see if one of them wants to trigger: 264 for (auto iModel : fApplicableModelList) { 265 if (iModel->AtRestModelTrigger(fFastTrack)) { 266 fFastStep.Initialize(fFastTrack); 267 fTriggedFastSimulationModel = iModel; 268 return true; 269 } 270 } 271 //-------------------------------------------- 272 // Nobody asks to gain control, returns false 273 //-------------------------------------------- 274 return false; 275 } 276 277 G4VParticleChange* G4FastSimulationManager::InvokeAtRestDoIt() 278 { 279 fTriggedFastSimulationModel->AtRestDoIt(fFastTrack, fFastStep); 280 return &fFastStep; 281 } 282 283 void G4FastSimulationManager::ListTitle() const 284 { 285 G4cout << fFastTrack.GetEnvelope()->GetName(); 286 if (fFastTrack.GetEnvelope()->GetWorldPhysical() 287 == G4TransportationManager::GetTransportationManager() 288 ->GetNavigatorForTracking() 289 ->GetWorldVolume()) 290 G4cout << " (mass geom.)"; 291 else 292 G4cout << " (// geom.)"; 293 } 294 295 void G4FastSimulationManager::ListModels() const 296 { 297 G4cout << "Current Models for the "; 298 ListTitle(); 299 G4cout << " envelope:\n"; 300 301 for (auto iModel : ModelList) 302 G4cout << " " << iModel->GetName() << "\n"; 303 304 for (auto iModel : fInactivatedModels) 305 G4cout << " " << iModel->GetName() << "(inactivated)\n"; 306 } 307 308 void G4FastSimulationManager::ListModels(const G4String& modelName) const 309 { 310 std::size_t iModel; 311 G4int titled = 0; 312 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable(); 313 314 // Active Models 315 for (iModel = 0; iModel < ModelList.size(); ++iModel) 316 if (ModelList[iModel]->GetName() == modelName || modelName == "all") { 317 if ((titled++) == 0) { 318 G4cout << "In the envelope "; 319 ListTitle(); 320 G4cout << ",\n"; 321 } 322 G4cout << " the model " << ModelList[iModel]->GetName() << " is applicable for :\n "; 323 324 G4int list_started = 0; 325 for (G4int iParticle = 0; iParticle < theParticleTable->entries(); iParticle++) 326 if (ModelList[iModel]->IsApplicable(*(theParticleTable->GetParticle(iParticle)))) { 327 if ((list_started++) != 0) G4cout << ", "; 328 G4cout << theParticleTable->GetParticle(iParticle)->GetParticleName(); 329 } 330 G4cout << G4endl; 331 } 332 333 // Inactive Models 334 for (iModel = 0; iModel < fInactivatedModels.size(); ++iModel) 335 if (fInactivatedModels[iModel]->GetName() == modelName || modelName == "all") { 336 if ((titled++) == 0) { 337 G4cout << "In the envelope "; 338 ListTitle(); 339 G4cout << ",\n"; 340 } 341 G4cout << " the model " << fInactivatedModels[iModel]->GetName() 342 << " (inactivated) is applicable for :\n "; 343 344 G4int list_started = 0; 345 for (G4int iParticle = 0; iParticle < theParticleTable->entries(); iParticle++) 346 if (fInactivatedModels[iModel]->IsApplicable(*(theParticleTable->GetParticle(iParticle)))) { 347 if ((list_started++) != 0) G4cout << ", "; 348 G4cout << theParticleTable->GetParticle(iParticle)->GetParticleName(); 349 } 350 G4cout << G4endl; 351 } 352 } 353 354 void G4FastSimulationManager::ListModels(const G4ParticleDefinition* particleDefinition) const 355 { 356 std::size_t iModel; 357 G4bool unique = true; 358 359 // Active Models 360 for (iModel = 0; iModel < ModelList.size(); ++iModel) 361 if (ModelList[iModel]->IsApplicable(*particleDefinition)) { 362 G4cout << "Envelope "; 363 ListTitle(); 364 G4cout << ", Model " << ModelList[iModel]->GetName() << "." << G4endl; 365 // -- Verify unicity of model attached to particleDefinition: 366 for (auto jModel = iModel + 1; jModel < ModelList.size(); jModel++) 367 if (ModelList[jModel]->IsApplicable(*particleDefinition)) unique = false; 368 } 369 370 // Inactive Models 371 for (iModel = 0; iModel < fInactivatedModels.size(); ++iModel) 372 if (fInactivatedModels[iModel]->IsApplicable(*particleDefinition)) { 373 G4cout << "Envelope "; 374 ListTitle(); 375 G4cout << ", Model " << fInactivatedModels[iModel]->GetName() << " (inactivated)." << G4endl; 376 } 377 378 if (!unique) { 379 G4ExceptionDescription ed; 380 ed << "Two or more active Models are available for the same particle type, in the same " 381 "envelope/region." 382 << G4endl; 383 G4Exception( 384 "G4FastSimulationManager::ListModels(const G4ParticleDefinition* particleDefinition) const", 385 "FastSim001", JustWarning, ed, "Models risk to exclude each other."); 386 } 387 unique = false; 388 } 389