Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 // 26 // 27 // ------------------------------------------- 27 // ------------------------------------------------------------------- 28 // 28 // 29 // GEANT4 Class file 29 // GEANT4 Class file 30 // 30 // 31 // 31 // 32 // File name: G4VMultipleScattering 32 // File name: G4VMultipleScattering 33 // 33 // 34 // Author: Vladimir Ivanchenko on base 34 // Author: Vladimir Ivanchenko on base of Laszlo Urban code 35 // 35 // 36 // Creation date: 25.03.2003 36 // Creation date: 25.03.2003 37 // 37 // 38 // Modifications: 38 // Modifications: 39 // 39 // 40 // 16-07-03 Use G4VMscModel interface (V.Ivanc 40 // 16-07-03 Use G4VMscModel interface (V.Ivanchenko) 41 // 03-11-03 Fix initialisation problem in Retr 41 // 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko) 42 // 04-11-03 Update PrintInfoDefinition (V.Ivan 42 // 04-11-03 Update PrintInfoDefinition (V.Ivanchenko) 43 // 01-03-04 SampleCosineTheta signature change 43 // 01-03-04 SampleCosineTheta signature changed 44 // 22-04-04 SampleCosineTheta signature change 44 // 22-04-04 SampleCosineTheta signature changed back to original 45 // 27-08-04 Add InitialiseForRun method (V.Iva 45 // 27-08-04 Add InitialiseForRun method (V.Ivanchneko) 46 // 08-11-04 Migration to new interface of Stor 46 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko) 47 // 11-03-05 Shift verbose level by 1 (V.Ivantc 47 // 11-03-05 Shift verbose level by 1 (V.Ivantchenko) 48 // 15-04-05 optimize internal interface (V.Iva 48 // 15-04-05 optimize internal interface (V.Ivanchenko) 49 // 15-04-05 remove boundary flag (V.Ivanchenko 49 // 15-04-05 remove boundary flag (V.Ivanchenko) 50 // 27-10-05 introduce virtual function MscStep 50 // 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko) 51 // 12-04-07 Add verbosity at destruction (V.Iv 51 // 12-04-07 Add verbosity at destruction (V.Ivanchenko) 52 // 27-10-07 Virtual functions moved to source 52 // 27-10-07 Virtual functions moved to source (V.Ivanchenko) 53 // 11-03-08 Set skin value does not effect ste 53 // 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko) 54 // 24-06-09 Removed hidden bin in G4PhysicsVec 54 // 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko) 55 // 04-06-13 Adoptation to MT mode (V.Ivanchenk 55 // 04-06-13 Adoptation to MT mode (V.Ivanchenko) 56 // 56 // 57 57 58 // ------------------------------------------- 58 // ------------------------------------------------------------------- 59 // 59 // 60 //....oooOO0OOooo........oooOO0OOooo........oo 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 61 //....oooOO0OOooo........oooOO0OOooo........oo 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 62 62 63 #include "G4VMultipleScattering.hh" 63 #include "G4VMultipleScattering.hh" 64 #include "G4PhysicalConstants.hh" 64 #include "G4PhysicalConstants.hh" 65 #include "G4SystemOfUnits.hh" 65 #include "G4SystemOfUnits.hh" 66 #include "G4LossTableManager.hh" 66 #include "G4LossTableManager.hh" 67 #include "G4MaterialCutsCouple.hh" 67 #include "G4MaterialCutsCouple.hh" 68 #include "G4Step.hh" 68 #include "G4Step.hh" 69 #include "G4ParticleDefinition.hh" 69 #include "G4ParticleDefinition.hh" 70 #include "G4VEmFluctuationModel.hh" 70 #include "G4VEmFluctuationModel.hh" 71 #include "G4UnitsTable.hh" 71 #include "G4UnitsTable.hh" 72 #include "G4ProductionCutsTable.hh" 72 #include "G4ProductionCutsTable.hh" 73 #include "G4Electron.hh" 73 #include "G4Electron.hh" 74 #include "G4GenericIon.hh" 74 #include "G4GenericIon.hh" 75 #include "G4TransportationManager.hh" 75 #include "G4TransportationManager.hh" 76 #include "G4SafetyHelper.hh" 76 #include "G4SafetyHelper.hh" 77 #include "G4ParticleTable.hh" 77 #include "G4ParticleTable.hh" 78 #include "G4ProcessVector.hh" 78 #include "G4ProcessVector.hh" 79 #include "G4ProcessManager.hh" 79 #include "G4ProcessManager.hh" 80 #include "G4LossTableBuilder.hh" 80 #include "G4LossTableBuilder.hh" 81 #include "G4EmTableUtil.hh" << 82 #include <iostream> 81 #include <iostream> 83 82 84 //....oooOO0OOooo........oooOO0OOooo........oo 83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 85 84 86 G4VMultipleScattering::G4VMultipleScattering(c 85 G4VMultipleScattering::G4VMultipleScattering(const G4String&, G4ProcessType) 87 : G4VContinuousDiscreteProcess("msc", fElect 86 : G4VContinuousDiscreteProcess("msc", fElectromagnetic), 88 fNewPosition(0.,0.,0.), 87 fNewPosition(0.,0.,0.), 89 fNewDirection(0.,0.,1.) 88 fNewDirection(0.,0.,1.) 90 { 89 { 91 theParameters = G4EmParameters::Instance(); 90 theParameters = G4EmParameters::Instance(); 92 SetVerboseLevel(1); 91 SetVerboseLevel(1); 93 SetProcessSubType(fMultipleScattering); 92 SetProcessSubType(fMultipleScattering); 94 93 95 lowestKinEnergy = 10*CLHEP::eV; 94 lowestKinEnergy = 10*CLHEP::eV; 96 95 97 geomMin = 0.05*CLHEP::nm; << 96 geomMin = 0.05*CLHEP::nm; 98 minDisplacement2 = geomMin*geomMin; 97 minDisplacement2 = geomMin*geomMin; 99 98 100 pParticleChange = &fParticleChange; 99 pParticleChange = &fParticleChange; 101 100 102 modelManager = new G4EmModelManager(); 101 modelManager = new G4EmModelManager(); 103 emManager = G4LossTableManager::Instance(); 102 emManager = G4LossTableManager::Instance(); 104 mscModels.reserve(2); 103 mscModels.reserve(2); 105 emManager->Register(this); 104 emManager->Register(this); 106 } 105 } 107 106 108 //....oooOO0OOooo........oooOO0OOooo........oo 107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 109 108 110 G4VMultipleScattering::~G4VMultipleScattering( 109 G4VMultipleScattering::~G4VMultipleScattering() 111 { 110 { 112 delete modelManager; 111 delete modelManager; 113 emManager->DeRegister(this); 112 emManager->DeRegister(this); 114 } 113 } 115 114 116 //....oooOO0OOooo........oooOO0OOooo........oo 115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 117 116 118 void G4VMultipleScattering::AddEmModel(G4int o << 117 void G4VMultipleScattering::AddEmModel(G4int order, G4VEmModel* ptr, 119 const G 118 const G4Region* region) 120 { 119 { 121 if(nullptr == ptr) { return; } 120 if(nullptr == ptr) { return; } 122 G4VEmFluctuationModel* fm = nullptr; 121 G4VEmFluctuationModel* fm = nullptr; 123 modelManager->AddEmModel(order, ptr, fm, reg 122 modelManager->AddEmModel(order, ptr, fm, region); 124 ptr->SetParticleChange(pParticleChange); 123 ptr->SetParticleChange(pParticleChange); 125 } 124 } 126 125 127 //....oooOO0OOooo........oooOO0OOooo........oo 126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 128 127 129 void G4VMultipleScattering::SetEmModel(G4VMscM 128 void G4VMultipleScattering::SetEmModel(G4VMscModel* ptr, G4int) 130 { 129 { 131 if(nullptr == ptr) { return; } 130 if(nullptr == ptr) { return; } 132 if(!mscModels.empty()) { 131 if(!mscModels.empty()) { 133 for(auto & msc : mscModels) { if(msc == pt 132 for(auto & msc : mscModels) { if(msc == ptr) { return; } } 134 } 133 } 135 mscModels.push_back(ptr); 134 mscModels.push_back(ptr); 136 } 135 } 137 136 138 //....oooOO0OOooo........oooOO0OOooo........oo 137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 139 138 140 void 139 void 141 G4VMultipleScattering::PreparePhysicsTable(con 140 G4VMultipleScattering::PreparePhysicsTable(const G4ParticleDefinition& part) 142 { 141 { >> 142 if(1 < verboseLevel) { >> 143 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for " >> 144 << GetProcessName() >> 145 << " and particle " << part.GetParticleName() >> 146 << G4endl; >> 147 } 143 G4bool master = emManager->IsMaster(); 148 G4bool master = emManager->IsMaster(); 144 if (nullptr == firstParticle) { firstParticl << 145 149 146 emManager->PreparePhysicsTable(&part, this); << 150 if(nullptr == firstParticle) { firstParticle = ∂ } >> 151 if(part.GetPDGMass() > CLHEP::GeV) { >> 152 // flag declears that mass scaling is applied >> 153 isIon = true; >> 154 } >> 155 >> 156 emManager->PreparePhysicsTable(&part, this, master); 147 currParticle = nullptr; 157 currParticle = nullptr; 148 158 >> 159 if(1 < verboseLevel) { >> 160 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for " >> 161 << GetProcessName() >> 162 << " and particle " << part.GetParticleName() >> 163 << " local particle " << firstParticle->GetParticleName() >> 164 << " isIon: " << isIon << " isMaster: " << master >> 165 << G4endl; >> 166 } >> 167 149 if(firstParticle == &part) { 168 if(firstParticle == &part) { 150 baseMat = emManager->GetTableBuilder()->Ge << 151 G4EmTableUtil::PrepareMscProcess(this, par << 152 stepLimit, facrange, << 153 latDisplacement, master, << 154 isIon, baseMat); << 155 169 >> 170 // initialise process >> 171 InitialiseProcess(firstParticle); >> 172 >> 173 // heavy particles >> 174 if(part.GetPDGMass() > CLHEP::MeV) { >> 175 stepLimit = theParameters->MscMuHadStepLimitType(); >> 176 facrange = theParameters->MscMuHadRangeFactor(); >> 177 latDisplacement = theParameters->MuHadLateralDisplacement(); >> 178 } else { >> 179 stepLimit = theParameters->MscStepLimitType(); >> 180 facrange = theParameters->MscRangeFactor(); >> 181 latDisplacement = theParameters->LateralDisplacement(); >> 182 } >> 183 if(master) { SetVerboseLevel(theParameters->Verbose()); } >> 184 else { SetVerboseLevel(theParameters->WorkerVerbose()); } >> 185 >> 186 // initialisation of models 156 numberOfModels = modelManager->NumberOfMod 187 numberOfModels = modelManager->NumberOfModels(); 157 currentModel = GetModelByIndex(0); << 188 /* >> 189 std::cout << "### G4VMultipleScattering::PreparePhysicsTable() for " >> 190 << GetProcessName() >> 191 << " and particle " << part.GetParticleName() >> 192 << " Nmodels= " << mscModels.size() << " " << this << std::endl; >> 193 */ >> 194 G4LossTableBuilder* bld = emManager->GetTableBuilder(); >> 195 baseMat = bld->GetBaseMaterialFlag(); >> 196 >> 197 for(G4int i=0; i<numberOfModels; ++i) { >> 198 G4VMscModel* msc = GetModelByIndex(i); >> 199 if(nullptr == msc) { continue; } >> 200 if(nullptr == currentModel) { currentModel = msc; } >> 201 msc->SetIonisation(nullptr, firstParticle); >> 202 msc->SetMasterThread(master); >> 203 msc->SetPolarAngleLimit(theParameters->MscThetaLimit()); >> 204 G4double emax = >> 205 std::min(msc->HighEnergyLimit(),theParameters->MaxKinEnergy()); >> 206 msc->SetHighEnergyLimit(emax); >> 207 msc->SetUseBaseMaterials(baseMat); >> 208 } >> 209 >> 210 modelManager->Initialise(firstParticle, G4Electron::Electron(), >> 211 1.0, verboseLevel); 158 212 159 if (nullptr == safetyHelper) { << 213 if(nullptr == safetyHelper) { 160 safetyHelper = G4TransportationManager:: 214 safetyHelper = G4TransportationManager::GetTransportationManager() 161 ->GetSafetyHelper(); << 215 ->GetSafetyHelper(); 162 safetyHelper->InitialiseHelper(); 216 safetyHelper->InitialiseHelper(); 163 } 217 } 164 } 218 } 165 } 219 } 166 220 167 //....oooOO0OOooo........oooOO0OOooo........oo 221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 168 222 169 void G4VMultipleScattering::BuildPhysicsTable( 223 void G4VMultipleScattering::BuildPhysicsTable(const G4ParticleDefinition& part) 170 { 224 { >> 225 const G4String& num = part.GetParticleName(); 171 G4bool master = emManager->IsMaster(); 226 G4bool master = emManager->IsMaster(); 172 << 227 if(1 < verboseLevel) { 173 if(firstParticle == &part) { << 228 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for " 174 emManager->BuildPhysicsTable(&part); << 229 << GetProcessName() >> 230 << " and particle " << num << " isIon: " << isIon >> 231 << " IsMaster: " << master << G4endl; 175 } 232 } 176 const G4VMultipleScattering* ptr = this; << 233 const G4VMultipleScattering* masterProcess = 177 if(!master) { << 234 static_cast<const G4VMultipleScattering*>(GetMasterProcess()); 178 ptr = static_cast<const G4VMultipleScatter << 235 >> 236 if(firstParticle == &part) { >> 237 /* >> 238 std::cout << "### G4VMultipleScattering::BuildPhysicsTable() for " >> 239 << GetProcessName() << " and particle " << num >> 240 << " IsMaster= " << G4LossTableManager::Instance()->IsMaster() >> 241 << " " << this << std::endl; >> 242 */ >> 243 emManager->BuildPhysicsTable(firstParticle); >> 244 >> 245 if(!master) { >> 246 // initialisation of models >> 247 /* >> 248 std::cout << "### G4VMultipleScattering::BuildPhysicsTable() for " >> 249 << GetProcessName() << " and particle " << num >> 250 << " Nmod= " << mscModels.size() << " NOT master" << std::endl; >> 251 */ >> 252 baseMat = masterProcess->UseBaseMaterial(); >> 253 for(G4int i=0; i<numberOfModels; ++i) { >> 254 G4VMscModel* msc = GetModelByIndex(i); >> 255 if(nullptr == msc) { continue; } >> 256 G4VMscModel* msc0 = masterProcess->GetModelByIndex(i); >> 257 msc->SetUseBaseMaterials(baseMat); >> 258 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false); >> 259 msc->InitialiseLocal(firstParticle, msc0); >> 260 } >> 261 } 179 } 262 } >> 263 // protection against double printout >> 264 if(theParameters->IsPrintLocked()) { return; } >> 265 >> 266 // explicitly defined printout by particle name >> 267 if(1 < verboseLevel || >> 268 (0 < verboseLevel && (num == "e-" || >> 269 num == "e+" || num == "mu+" || >> 270 num == "mu-" || num == "proton"|| >> 271 num == "pi+" || num == "pi-" || >> 272 num == "kaon+" || num == "kaon-" || >> 273 num == "alpha" || num == "anti_proton" || >> 274 num == "GenericIon" || num == "alpha+" || >> 275 num == "alpha++" ))) >> 276 { >> 277 StreamInfo(G4cout, part); >> 278 } 180 279 181 G4EmTableUtil::BuildMscProcess(this, ptr, pa << 280 if(1 < verboseLevel) { 182 numberOfModels, master); << 281 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for " >> 282 << GetProcessName() >> 283 << " and particle " << num << G4endl; >> 284 } 183 } 285 } 184 286 185 //....oooOO0OOooo........oooOO0OOooo........oo 287 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 186 288 187 void G4VMultipleScattering::StreamInfo(std::os 289 void G4VMultipleScattering::StreamInfo(std::ostream& outFile, 188 const G4ParticleDefinition& 290 const G4ParticleDefinition& part, G4bool rst) const 189 { 291 { 190 G4String indent = (rst ? " " : ""); 292 G4String indent = (rst ? " " : ""); 191 outFile << G4endl << indent << GetProcessNam 293 outFile << G4endl << indent << GetProcessName() << ": "; 192 if (!rst) outFile << " for " << part.GetPart 294 if (!rst) outFile << " for " << part.GetParticleName(); 193 outFile << " SubType= " << GetProcessSubTy 295 outFile << " SubType= " << GetProcessSubType() << G4endl; 194 modelManager->DumpModelList(outFile, verbose 296 modelManager->DumpModelList(outFile, verboseLevel); 195 } 297 } 196 298 197 //....oooOO0OOooo........oooOO0OOooo........oo 299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 198 300 199 void G4VMultipleScattering::StartTracking(G4Tr 301 void G4VMultipleScattering::StartTracking(G4Track* track) 200 { 302 { 201 G4VEnergyLossProcess* eloss = nullptr; 303 G4VEnergyLossProcess* eloss = nullptr; 202 if(track->GetParticleDefinition() != currPar 304 if(track->GetParticleDefinition() != currParticle) { 203 currParticle = track->GetParticleDefinitio 305 currParticle = track->GetParticleDefinition(); 204 fIonisation = emManager->GetEnergyLossProc 306 fIonisation = emManager->GetEnergyLossProcess(currParticle); 205 eloss = fIonisation; 307 eloss = fIonisation; 206 } 308 } >> 309 /* >> 310 G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels >> 311 << " " << currParticle->GetParticleName() >> 312 << " E(MeV)= " << track->GetKineticEnergy() >> 313 << " Ion= " << fIonisation << " IsMaster= " >> 314 << G4LossTableManager::Instance()->IsMaster() >> 315 << G4endl; >> 316 */ 207 for(G4int i=0; i<numberOfModels; ++i) { 317 for(G4int i=0; i<numberOfModels; ++i) { 208 G4VMscModel* msc = GetModelByIndex(i); 318 G4VMscModel* msc = GetModelByIndex(i); >> 319 /* >> 320 G4cout << "Next model " << msc >> 321 << " Emin= " << msc->LowEnergyLimit() >> 322 << " Emax= " << msc->HighEnergyLimit() >> 323 << " Eact= " << msc->LowEnergyActivationLimit() << G4endl; >> 324 */ 209 msc->StartTracking(track); 325 msc->StartTracking(track); 210 if(nullptr != eloss) { 326 if(nullptr != eloss) { 211 msc->SetIonisation(eloss, currParticle); 327 msc->SetIonisation(eloss, currParticle); 212 } 328 } 213 } 329 } 214 } 330 } 215 331 216 //....oooOO0OOooo........oooOO0OOooo........oo 332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 217 333 218 G4double G4VMultipleScattering::AlongStepGetPh 334 G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength( 219 const G4Track& tr 335 const G4Track& track, 220 G4double, 336 G4double, 221 G4double currentM 337 G4double currentMinimalStep, 222 G4double&, 338 G4double&, 223 G4GPILSelection* 339 G4GPILSelection* selection) 224 { 340 { 225 // get Step limit proposed by the process 341 // get Step limit proposed by the process 226 *selection = NotCandidateForSelection; 342 *selection = NotCandidateForSelection; 227 physStepLimit = gPathLength = tPathLength = 343 physStepLimit = gPathLength = tPathLength = currentMinimalStep; 228 344 229 G4double ekin = track.GetKineticEnergy(); 345 G4double ekin = track.GetKineticEnergy(); 230 /* 346 /* 231 G4cout << "MSC::AlongStepGPIL: Ekin= " << ek 347 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin 232 << " " << currParticle->GetParticleN 348 << " " << currParticle->GetParticleName() 233 << " currMod " << currentModel 349 << " currMod " << currentModel 234 << G4endl; 350 << G4endl; 235 */ 351 */ 236 // isIon flag is used only to select a model 352 // isIon flag is used only to select a model 237 if(isIon) { 353 if(isIon) { 238 ekin *= proton_mass_c2/track.GetParticleDe 354 ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass(); 239 } 355 } 240 const G4MaterialCutsCouple* couple = track.G 356 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple(); 241 357 242 // select new model, static cast is possible << 358 // select new model 243 if(1 < numberOfModels) { 359 if(1 < numberOfModels) { 244 currentModel = << 360 currentModel = static_cast<G4VMscModel*>(SelectModel(ekin,couple->GetIndex())); 245 static_cast<G4VMscModel*>(SelectModel(ek << 246 } 361 } 247 currentModel->SetCurrentCouple(couple); 362 currentModel->SetCurrentCouple(couple); 248 // msc is active is model is active, energy 363 // msc is active is model is active, energy above the limit, 249 // and step size is above the limit; 364 // and step size is above the limit; 250 // if it is active msc may limit the step 365 // if it is active msc may limit the step 251 if(currentModel->IsActive(ekin) && tPathLeng 366 if(currentModel->IsActive(ekin) && tPathLength > geomMin 252 && ekin >= lowestKinEnergy) { 367 && ekin >= lowestKinEnergy) { 253 isActive = true; 368 isActive = true; 254 tPathLength = 369 tPathLength = 255 currentModel->ComputeTruePathLengthLimit 370 currentModel->ComputeTruePathLengthLimit(track, gPathLength); 256 if (tPathLength < physStepLimit) { 371 if (tPathLength < physStepLimit) { 257 *selection = CandidateForSelection; 372 *selection = CandidateForSelection; 258 } 373 } 259 } else { << 374 } else { isActive = false; } 260 isActive = false; << 261 gPathLength = DBL_MAX; << 262 } << 263 375 264 //if(currParticle->GetPDGMass() > GeV) 376 //if(currParticle->GetPDGMass() > GeV) 265 /* 377 /* 266 G4cout << "MSC::AlongStepGPIL: Ekin= " << ek 378 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin 267 << " " << currParticle->GetParticleN 379 << " " << currParticle->GetParticleName() 268 << " gPathLength= " << gPathLength 380 << " gPathLength= " << gPathLength 269 << " tPathLength= " << tPathLength 381 << " tPathLength= " << tPathLength 270 << " currentMinimalStep= " << current 382 << " currentMinimalStep= " << currentMinimalStep 271 << " isActive " << isActive << G4endl 383 << " isActive " << isActive << G4endl; 272 */ 384 */ 273 return gPathLength; 385 return gPathLength; 274 } 386 } 275 387 276 //....oooOO0OOooo........oooOO0OOooo........oo 388 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 277 389 278 G4double 390 G4double 279 G4VMultipleScattering::PostStepGetPhysicalInte 391 G4VMultipleScattering::PostStepGetPhysicalInteractionLength( 280 const G4Track&, G4double, G4Forc 392 const G4Track&, G4double, G4ForceCondition* condition) 281 { 393 { 282 *condition = NotForced; 394 *condition = NotForced; 283 return DBL_MAX; 395 return DBL_MAX; 284 } 396 } 285 397 286 //....oooOO0OOooo........oooOO0OOooo........oo 398 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 287 399 288 G4VParticleChange* 400 G4VParticleChange* 289 G4VMultipleScattering::AlongStepDoIt(const G4T 401 G4VMultipleScattering::AlongStepDoIt(const G4Track& track, const G4Step& step) 290 { 402 { 291 fParticleChange.InitialiseMSC(track, step); << 403 fParticleChange.ProposeMomentumDirection( 292 fNewPosition = fParticleChange.GetProposedPo << 404 step.GetPostStepPoint()->GetMomentumDirection()); >> 405 fNewPosition = step.GetPostStepPoint()->GetPosition(); >> 406 fParticleChange.ProposePosition(fNewPosition); 293 fPositionChanged = false; 407 fPositionChanged = false; 294 408 295 G4double geomLength = step.GetStepLength(); 409 G4double geomLength = step.GetStepLength(); 296 410 297 // very small step - no msc 411 // very small step - no msc 298 if(!isActive) { 412 if(!isActive) { 299 tPathLength = geomLength; 413 tPathLength = geomLength; 300 414 301 // sample msc 415 // sample msc 302 } else { 416 } else { 303 G4double range = 417 G4double range = 304 currentModel->GetRange(currParticle,trac 418 currentModel->GetRange(currParticle,track.GetKineticEnergy(), 305 track.GetMaterial 419 track.GetMaterialCutsCouple()); 306 420 307 tPathLength = currentModel->ComputeTrueSte 421 tPathLength = currentModel->ComputeTrueStepLength(geomLength); 308 422 309 /* 423 /* 310 if(currParticle->GetPDGMass() > 0.9*GeV) 424 if(currParticle->GetPDGMass() > 0.9*GeV) 311 G4cout << "G4VMsc::AlongStepDoIt: GeomLeng 425 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= " 312 << geomLength 426 << geomLength 313 << " tPathLength= " << tPathLength 427 << " tPathLength= " << tPathLength 314 << " physStepLimit= " << physStepLi 428 << " physStepLimit= " << physStepLimit 315 << " dr= " << range - tPathLength 429 << " dr= " << range - tPathLength 316 << " ekin= " << track.GetKineticEne 430 << " ekin= " << track.GetKineticEnergy() << G4endl; 317 */ 431 */ 318 // protection against wrong t->g->t conver 432 // protection against wrong t->g->t conversion 319 tPathLength = std::min(tPathLength, physSt 433 tPathLength = std::min(tPathLength, physStepLimit); 320 434 321 // do not sample scattering at the last or 435 // do not sample scattering at the last or at a small step 322 if(tPathLength < range && tPathLength > ge 436 if(tPathLength < range && tPathLength > geomMin) { 323 437 324 static const G4double minSafety = 1.20*C 438 static const G4double minSafety = 1.20*CLHEP::nm; 325 static const G4double sFact = 0.99; 439 static const G4double sFact = 0.99; 326 440 327 G4ThreeVector displacement = currentMode 441 G4ThreeVector displacement = currentModel->SampleScattering( 328 step.GetPostStepPoint()->GetMomentumDi 442 step.GetPostStepPoint()->GetMomentumDirection(),minSafety); 329 443 330 G4double r2 = displacement.mag2(); 444 G4double r2 = displacement.mag2(); 331 //G4cout << " R= " << sqrt(r2) << " R 445 //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2) 332 // << " flag= " << fDispBeyondSafety 446 // << " flag= " << fDispBeyondSafety << G4endl; 333 if(r2 > minDisplacement2) { 447 if(r2 > minDisplacement2) { 334 448 335 fPositionChanged = true; 449 fPositionChanged = true; 336 G4double dispR = std::sqrt(r2); 450 G4double dispR = std::sqrt(r2); 337 G4double postSafety = 451 G4double postSafety = 338 sFact*safetyHelper->ComputeSafety(fN 452 sFact*safetyHelper->ComputeSafety(fNewPosition, dispR); 339 //G4cout<<" R= "<< dispR<<" postSaf 453 //G4cout<<" R= "<< dispR<<" postSafety= "<<postSafety<<G4endl; 340 454 341 // far away from geometry boundary 455 // far away from geometry boundary 342 if(postSafety > 0.0 && dispR <= postSa 456 if(postSafety > 0.0 && dispR <= postSafety) { 343 fNewPosition += displacement; 457 fNewPosition += displacement; 344 458 345 //near the boundary 459 //near the boundary 346 } else { 460 } else { 347 // displaced point is definitely wit 461 // displaced point is definitely within the volume 348 //G4cout<<" R= "<<dispR<<" postSa 462 //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl; 349 if(dispR < postSafety) { 463 if(dispR < postSafety) { 350 fNewPosition += displacement; 464 fNewPosition += displacement; 351 465 352 // reduced displacement 466 // reduced displacement 353 } else if(postSafety > geomMin) { 467 } else if(postSafety > geomMin) { 354 fNewPosition += displacement*(post 468 fNewPosition += displacement*(postSafety/dispR); 355 469 356 // very small postSafety 470 // very small postSafety 357 } else { 471 } else { 358 fPositionChanged = false; 472 fPositionChanged = false; 359 } 473 } 360 } 474 } 361 if(fPositionChanged) { 475 if(fPositionChanged) { 362 safetyHelper->ReLocateWithinVolume(f 476 safetyHelper->ReLocateWithinVolume(fNewPosition); 363 fParticleChange.ProposePosition(fNew 477 fParticleChange.ProposePosition(fNewPosition); 364 } 478 } 365 } 479 } 366 } 480 } 367 } 481 } 368 fParticleChange.ProposeTrueStepLength(tPathL 482 fParticleChange.ProposeTrueStepLength(tPathLength); 369 return &fParticleChange; 483 return &fParticleChange; 370 } 484 } 371 485 372 //....oooOO0OOooo........oooOO0OOooo........oo 486 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 373 487 >> 488 G4VParticleChange* >> 489 G4VMultipleScattering::PostStepDoIt(const G4Track& track, const G4Step&) >> 490 { >> 491 fParticleChange.Initialize(track); >> 492 return &fParticleChange; >> 493 } >> 494 >> 495 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... >> 496 374 G4double G4VMultipleScattering::GetContinuousS 497 G4double G4VMultipleScattering::GetContinuousStepLimit( 375 const G 498 const G4Track& track, 376 G4doubl 499 G4double previousStepSize, 377 G4doubl 500 G4double currentMinimalStep, 378 G4doubl 501 G4double& currentSafety) 379 { 502 { 380 G4GPILSelection selection = NotCandidateForS 503 G4GPILSelection selection = NotCandidateForSelection; 381 G4double x = AlongStepGetPhysicalInteraction 504 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize, 382 505 currentMinimalStep, 383 506 currentSafety, 384 507 &selection); 385 return x; 508 return x; 386 } 509 } 387 510 388 //....oooOO0OOooo........oooOO0OOooo........oo 511 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 389 512 390 G4double G4VMultipleScattering::ContinuousStep 513 G4double G4VMultipleScattering::ContinuousStepLimit( 391 const G 514 const G4Track& track, 392 G4doubl 515 G4double previousStepSize, 393 G4doubl 516 G4double currentMinimalStep, 394 G4doubl 517 G4double& currentSafety) 395 { 518 { 396 return GetContinuousStepLimit(track,previous 519 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep, 397 currentSafety) 520 currentSafety); 398 } 521 } 399 522 400 //....oooOO0OOooo........oooOO0OOooo........oo 523 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 401 524 402 G4double G4VMultipleScattering::GetMeanFreePat 525 G4double G4VMultipleScattering::GetMeanFreePath( 403 const G4Track&, G4double, G4Forc 526 const G4Track&, G4double, G4ForceCondition* condition) 404 { 527 { 405 *condition = Forced; 528 *condition = Forced; 406 return DBL_MAX; 529 return DBL_MAX; 407 } 530 } 408 531 409 //....oooOO0OOooo........oooOO0OOooo........oo 532 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 410 533 411 G4bool 534 G4bool 412 G4VMultipleScattering::StorePhysicsTable(const 535 G4VMultipleScattering::StorePhysicsTable(const G4ParticleDefinition* part, 413 const 536 const G4String& directory, 414 G4boo 537 G4bool ascii) 415 { 538 { 416 G4bool yes = true; 539 G4bool yes = true; 417 if(part != firstParticle || !emManager->IsMa << 540 if(part != firstParticle) { return yes; } 418 << 541 const G4VMultipleScattering* masterProcess = 419 return G4EmTableUtil::StoreMscTable(this, pa << 542 static_cast<const G4VMultipleScattering*>(GetMasterProcess()); 420 numberOfModels, verboseLevel, << 543 if(nullptr != masterProcess && masterProcess != this) { return yes; } 421 ascii); << 544 >> 545 G4int nmod = modelManager->NumberOfModels(); >> 546 static const G4String ss[4] = {"1","2","3","4"}; >> 547 for(G4int i=0; i<nmod; ++i) { >> 548 G4VEmModel* msc = modelManager->GetModel(i); >> 549 if(nullptr == msc) { continue; } >> 550 yes = true; >> 551 G4PhysicsTable* table = msc->GetCrossSectionTable(); >> 552 if (nullptr != table) { >> 553 G4int j = std::min(i,3); >> 554 G4String name = >> 555 GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii); >> 556 yes = table->StorePhysicsTable(name,ascii); >> 557 >> 558 if ( yes ) { >> 559 if ( verboseLevel>0 ) { >> 560 G4cout << "Physics table are stored for " >> 561 << part->GetParticleName() >> 562 << " and process " << GetProcessName() >> 563 << " with a name <" << name << "> " << G4endl; >> 564 } >> 565 } else { >> 566 G4cout << "Fail to store Physics Table for " >> 567 << part->GetParticleName() >> 568 << " and process " << GetProcessName() >> 569 << " in the directory <" << directory >> 570 << "> " << G4endl; >> 571 } >> 572 } >> 573 } >> 574 return yes; 422 } 575 } 423 576 424 //....oooOO0OOooo........oooOO0OOooo........oo 577 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 425 578 426 G4bool 579 G4bool 427 G4VMultipleScattering::RetrievePhysicsTable(co 580 G4VMultipleScattering::RetrievePhysicsTable(const G4ParticleDefinition*, 428 co 581 const G4String&, 429 G4 582 G4bool) 430 { 583 { 431 return true; 584 return true; 432 } 585 } 433 586 434 //....oooOO0OOooo........oooOO0OOooo........oo 587 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 435 588 >> 589 void G4VMultipleScattering::SetIonisation(G4VEnergyLossProcess* p) >> 590 { >> 591 for(auto & msc : mscModels) { >> 592 if(nullptr != msc) { msc->SetIonisation(p, firstParticle); } >> 593 } >> 594 } >> 595 >> 596 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... >> 597 436 void G4VMultipleScattering::ProcessDescription 598 void G4VMultipleScattering::ProcessDescription(std::ostream& outFile) const 437 { 599 { 438 if(nullptr != firstParticle) { << 600 if(firstParticle) { 439 StreamInfo(outFile, *firstParticle, true); 601 StreamInfo(outFile, *firstParticle, true); 440 } 602 } 441 } 603 } 442 604 443 //....oooOO0OOooo........oooOO0OOooo........oo 605 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 444 606 445 607