Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/utils/src/G4VMultipleScattering.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // GEANT4 Class file
 30 //
 31 //
 32 // File name:     G4VMultipleScattering
 33 //
 34 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
 35 //
 36 // Creation date: 25.03.2003
 37 //
 38 // Modifications:
 39 //
 40 // 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
 41 // 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
 42 // 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
 43 // 01-03-04 SampleCosineTheta signature changed
 44 // 22-04-04 SampleCosineTheta signature changed back to original
 45 // 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
 46 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
 47 // 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
 48 // 15-04-05 optimize internal interface (V.Ivanchenko)
 49 // 15-04-05 remove boundary flag (V.Ivanchenko)
 50 // 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
 51 // 12-04-07 Add verbosity at destruction (V.Ivanchenko)
 52 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
 53 // 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
 54 // 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
 55 // 04-06-13 Adoptation to MT mode (V.Ivanchenko)
 56 //
 57 
 58 // -------------------------------------------------------------------
 59 //
 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 62 
 63 #include "G4VMultipleScattering.hh"
 64 #include "G4PhysicalConstants.hh"
 65 #include "G4SystemOfUnits.hh"
 66 #include "G4LossTableManager.hh"
 67 #include "G4MaterialCutsCouple.hh"
 68 #include "G4Step.hh"
 69 #include "G4ParticleDefinition.hh"
 70 #include "G4VEmFluctuationModel.hh"
 71 #include "G4UnitsTable.hh"
 72 #include "G4ProductionCutsTable.hh"
 73 #include "G4Electron.hh"
 74 #include "G4GenericIon.hh"
 75 #include "G4TransportationManager.hh"
 76 #include "G4SafetyHelper.hh"
 77 #include "G4ParticleTable.hh"
 78 #include "G4ProcessVector.hh"
 79 #include "G4ProcessManager.hh"
 80 #include "G4LossTableBuilder.hh"
 81 #include "G4EmTableUtil.hh"
 82 #include <iostream>
 83 
 84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 85 
 86 G4VMultipleScattering::G4VMultipleScattering(const G4String&, G4ProcessType)
 87   : G4VContinuousDiscreteProcess("msc", fElectromagnetic),
 88   fNewPosition(0.,0.,0.),
 89   fNewDirection(0.,0.,1.)
 90 {
 91   theParameters = G4EmParameters::Instance();
 92   SetVerboseLevel(1);
 93   SetProcessSubType(fMultipleScattering);
 94 
 95   lowestKinEnergy = 10*CLHEP::eV;
 96 
 97   geomMin = 0.05*CLHEP::nm;
 98   minDisplacement2 = geomMin*geomMin;
 99 
100   pParticleChange = &fParticleChange;
101 
102   modelManager = new G4EmModelManager();
103   emManager = G4LossTableManager::Instance();
104   mscModels.reserve(2);
105   emManager->Register(this);
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109 
110 G4VMultipleScattering::~G4VMultipleScattering()
111 {
112   delete modelManager;
113   emManager->DeRegister(this);
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117 
118 void G4VMultipleScattering::AddEmModel(G4int order, G4VMscModel* ptr,
119                                        const G4Region* region)
120 {
121   if(nullptr == ptr) { return; }
122   G4VEmFluctuationModel* fm = nullptr;
123   modelManager->AddEmModel(order, ptr, fm, region);
124   ptr->SetParticleChange(pParticleChange);
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128 
129 void G4VMultipleScattering::SetEmModel(G4VMscModel* ptr, G4int)
130 {
131   if(nullptr == ptr) { return; }
132   if(!mscModels.empty()) { 
133     for(auto & msc : mscModels) { if(msc == ptr) { return; } } 
134   }
135   mscModels.push_back(ptr);
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
140 void 
141 G4VMultipleScattering::PreparePhysicsTable(const G4ParticleDefinition& part)
142 {
143   G4bool master = emManager->IsMaster();
144   if (nullptr == firstParticle) { firstParticle = &part; }
145 
146   emManager->PreparePhysicsTable(&part, this);
147   currParticle = nullptr;
148 
149   if(firstParticle == &part) {
150     baseMat = emManager->GetTableBuilder()->GetBaseMaterialFlag();
151     G4EmTableUtil::PrepareMscProcess(this, part, modelManager,
152              stepLimit, facrange,
153              latDisplacement, master,
154              isIon, baseMat);
155 
156     numberOfModels = modelManager->NumberOfModels();
157     currentModel = GetModelByIndex(0);
158 
159     if (nullptr == safetyHelper) {
160       safetyHelper = G4TransportationManager::GetTransportationManager()
161   ->GetSafetyHelper();
162       safetyHelper->InitialiseHelper();
163     }
164   }
165 }
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168 
169 void G4VMultipleScattering::BuildPhysicsTable(const G4ParticleDefinition& part)
170 {
171   G4bool master = emManager->IsMaster();
172 
173   if(firstParticle == &part) {
174     emManager->BuildPhysicsTable(&part);
175   }
176   const G4VMultipleScattering* ptr = this;
177   if(!master) {
178     ptr = static_cast<const G4VMultipleScattering*>(GetMasterProcess());
179   }
180 
181   G4EmTableUtil::BuildMscProcess(this, ptr, part, firstParticle,
182          numberOfModels, master);
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186 
187 void G4VMultipleScattering::StreamInfo(std::ostream& outFile, 
188                   const G4ParticleDefinition& part, G4bool rst) const
189 {
190   G4String indent = (rst ? "  " : "");
191   outFile << G4endl << indent << GetProcessName() << ": ";
192   if (!rst) outFile << " for " << part.GetParticleName();
193   outFile  << "  SubType= " << GetProcessSubType() << G4endl;
194   modelManager->DumpModelList(outFile, verboseLevel);
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198 
199 void G4VMultipleScattering::StartTracking(G4Track* track)
200 {
201   G4VEnergyLossProcess* eloss = nullptr;
202   if(track->GetParticleDefinition() != currParticle) {
203     currParticle = track->GetParticleDefinition();
204     fIonisation = emManager->GetEnergyLossProcess(currParticle);
205     eloss = fIonisation;
206   }
207   for(G4int i=0; i<numberOfModels; ++i) {
208     G4VMscModel* msc = GetModelByIndex(i);
209     msc->StartTracking(track);
210     if(nullptr != eloss) {
211       msc->SetIonisation(eloss, currParticle);
212     }
213   }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
218 G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(
219                              const G4Track& track,
220                              G4double,
221                              G4double currentMinimalStep,
222                              G4double&,
223                              G4GPILSelection* selection)
224 {
225   // get Step limit proposed by the process
226   *selection = NotCandidateForSelection;
227   physStepLimit = gPathLength = tPathLength = currentMinimalStep;
228 
229   G4double ekin = track.GetKineticEnergy();
230   /*
231   G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
232          << "  " << currParticle->GetParticleName() 
233          << " currMod " << currentModel 
234          << G4endl;
235   */
236   // isIon flag is used only to select a model
237   if(isIon) { 
238     ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass(); 
239   }
240   const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
241   
242   // select new model, static cast is possible in this class
243   if(1 < numberOfModels) {
244     currentModel =
245       static_cast<G4VMscModel*>(SelectModel(ekin,couple->GetIndex()));
246   }
247   currentModel->SetCurrentCouple(couple);
248   // msc is active is model is active, energy above the limit,
249   // and step size is above the limit;
250   // if it is active msc may limit the step
251   if(currentModel->IsActive(ekin) && tPathLength > geomMin
252      && ekin >= lowestKinEnergy) {
253     isActive = true;
254     tPathLength = 
255       currentModel->ComputeTruePathLengthLimit(track, gPathLength);
256     if (tPathLength < physStepLimit) { 
257       *selection = CandidateForSelection; 
258     }
259   } else { 
260     isActive = false;
261     gPathLength = DBL_MAX;
262   }
263   
264   //if(currParticle->GetPDGMass() > GeV)    
265   /*
266   G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
267          << "  " << currParticle->GetParticleName()
268          << " gPathLength= " << gPathLength
269          << " tPathLength= " << tPathLength
270          << " currentMinimalStep= " << currentMinimalStep
271          << " isActive " << isActive << G4endl;
272   */
273   return gPathLength;
274 }
275 
276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
277 
278 G4double 
279 G4VMultipleScattering::PostStepGetPhysicalInteractionLength(
280               const G4Track&, G4double, G4ForceCondition* condition)
281 {
282   *condition = NotForced;
283   return DBL_MAX;
284 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287 
288 G4VParticleChange* 
289 G4VMultipleScattering::AlongStepDoIt(const G4Track& track, const G4Step& step)
290 {
291   fParticleChange.InitialiseMSC(track, step);
292   fNewPosition = fParticleChange.GetProposedPosition();
293   fPositionChanged = false;
294 
295   G4double geomLength = step.GetStepLength();
296 
297   // very small step - no msc
298   if(!isActive) {
299     tPathLength = geomLength;
300 
301     // sample msc
302   } else {
303     G4double range = 
304       currentModel->GetRange(currParticle,track.GetKineticEnergy(),
305                              track.GetMaterialCutsCouple());
306 
307     tPathLength = currentModel->ComputeTrueStepLength(geomLength);
308   
309     /*    
310     if(currParticle->GetPDGMass() > 0.9*GeV)    
311     G4cout << "G4VMsc::AlongStepDoIt: GeomLength= " 
312            << geomLength 
313            << " tPathLength= " << tPathLength
314            << " physStepLimit= " << physStepLimit
315            << " dr= " << range - tPathLength
316            << " ekin= " << track.GetKineticEnergy() << G4endl;
317     */
318     // protection against wrong t->g->t conversion
319     tPathLength = std::min(tPathLength, physStepLimit);
320 
321     // do not sample scattering at the last or at a small step
322     if(tPathLength < range && tPathLength > geomMin) {
323 
324       static const G4double minSafety = 1.20*CLHEP::nm;
325       static const G4double sFact = 0.99;
326 
327       G4ThreeVector displacement = currentModel->SampleScattering(
328         step.GetPostStepPoint()->GetMomentumDirection(),minSafety);
329 
330       G4double r2 = displacement.mag2();
331       //G4cout << "    R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
332       //     << " flag= " << fDispBeyondSafety << G4endl;
333       if(r2 > minDisplacement2) {
334 
335         fPositionChanged = true;
336         G4double dispR = std::sqrt(r2);
337         G4double postSafety = 
338           sFact*safetyHelper->ComputeSafety(fNewPosition, dispR); 
339         //G4cout<<"    R= "<< dispR<<" postSafety= "<<postSafety<<G4endl;
340 
341         // far away from geometry boundary
342         if(postSafety > 0.0 && dispR <= postSafety) {
343           fNewPosition += displacement;
344 
345           //near the boundary
346         } else {
347           // displaced point is definitely within the volume
348           //G4cout<<"    R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
349           if(dispR < postSafety) {
350             fNewPosition += displacement;
351 
352             // reduced displacement
353           } else if(postSafety > geomMin) {
354             fNewPosition += displacement*(postSafety/dispR); 
355 
356             // very small postSafety
357           } else {
358             fPositionChanged = false;
359           }
360         }
361         if(fPositionChanged) { 
362           safetyHelper->ReLocateWithinVolume(fNewPosition);
363           fParticleChange.ProposePosition(fNewPosition); 
364         }
365       }
366     }
367   }
368   fParticleChange.ProposeTrueStepLength(tPathLength);
369   return &fParticleChange;
370 }
371 
372 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
373 
374 G4double G4VMultipleScattering::GetContinuousStepLimit(
375                                        const G4Track& track,
376                                        G4double previousStepSize,
377                                        G4double currentMinimalStep,
378                                        G4double& currentSafety)
379 {
380   G4GPILSelection selection = NotCandidateForSelection;
381   G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
382                                                      currentMinimalStep,
383                                                      currentSafety, 
384                                                      &selection);
385   return x;
386 }
387 
388 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
389 
390 G4double G4VMultipleScattering::ContinuousStepLimit(
391                                        const G4Track& track,
392                                        G4double previousStepSize,
393                                        G4double currentMinimalStep,
394                                        G4double& currentSafety)
395 {
396   return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
397                                 currentSafety);
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401 
402 G4double G4VMultipleScattering::GetMeanFreePath(
403               const G4Track&, G4double, G4ForceCondition* condition)
404 {
405   *condition = Forced;
406   return DBL_MAX;
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
410 
411 G4bool 
412 G4VMultipleScattering::StorePhysicsTable(const G4ParticleDefinition* part,
413                                          const G4String& directory,
414                                          G4bool ascii)
415 {
416   G4bool yes = true;
417   if(part != firstParticle || !emManager->IsMaster()) { return yes; }
418 
419   return G4EmTableUtil::StoreMscTable(this, part, directory,
420               numberOfModels, verboseLevel,
421                                       ascii);
422 }
423 
424 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
425 
426 G4bool 
427 G4VMultipleScattering::RetrievePhysicsTable(const G4ParticleDefinition*,
428                                             const G4String&,
429                                             G4bool)
430 {
431   return true;
432 }
433 
434 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
435 
436 void G4VMultipleScattering::ProcessDescription(std::ostream& outFile) const
437 {
438   if(nullptr != firstParticle) {
439     StreamInfo(outFile, *firstParticle, true);
440   }
441 }
442 
443 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
444 
445