Geant4-11
G4DNABrownianTransportation.cc
Go to the documentation of this file.
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// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
28//
29// WARNING : This class is released as a prototype.
30// It might strongly evolve or even disapear in the next releases.
31//
32// History:
33// -----------
34// 10 Oct 2011 M.Karamitros created
35//
36// -------------------------------------------------------------------
37
40
41#include <CLHEP/Random/Stat.h>
42
44
45#include <G4Scheduler.hh>
47#include "G4SystemOfUnits.hh"
48#include "Randomize.hh"
49#include "G4Molecule.hh"
50#include "G4RandomDirection.hh"
51#include "G4ParticleTable.hh"
52#include "G4SafetyHelper.hh"
54#include "G4UnitsTable.hh"
55#include "G4NistManager.hh"
57#include "G4ITNavigator.hh"
58#include "G4ITSafetyHelper.hh" // Not used yet
61
62using namespace std;
63
64#ifndef State
65#define State(theXInfo) (GetState<G4ITBrownianState>()->theXInfo)
66#endif
67
68//#ifndef State
69//#define State(theXInfo) (fTransportationState->theXInfo)
70//#endif
71
72//#define USE_COLOR 1
73
74#ifdef USE_COLOR
75#define RED "\033[0;31m"
76#define LIGHT_RED "\33[1;31m"
77#define GREEN "\033[32;40m"
78#define GREEN_ON_BLUE "\033[1;32;44m"
79#define RESET_COLOR "\033[0m"
80#else
81#define RED ""
82#define LIGHT_RED ""
83#define GREEN ""
84#define GREEN_ON_BLUE ""
85#define RESET_COLOR ""
86#endif
87
88//#define DEBUG_MEM 1
89
90#ifdef DEBUG_MEM
91#include "G4MemStat.hh"
92using namespace G4MemStat;
94#endif
95
97{
99}
100
102{
103 return CLHEP::HepStat::inverseErf(1. - x);
104}
105
106//static double Erf(double x)
107//{
108// return CLHEP::HepStat::erf(x);
109//}
110
112{
113 return 1 - CLHEP::HepStat::erf(1. - x);
114}
115
116#ifndef State
117#define State(theXInfo) (GetState<G4ITTransportationState>()->theXInfo)
118#endif
119
121 G4int verbosity) :
122 G4ITTransportation(aName, verbosity)
123{
124
125 fVerboseLevel = 0;
126
127 fpState.reset(new G4ITBrownianState());
128
129 //ctor
131
133 fpWaterDensity = 0;
134
137 fSpeedMeUp = true;
138
141}
142
144{
146}
147
149 G4ITTransportation(right)
150{
151 //copy ctor
156 fNistWater = right.fNistWater;
159 fSpeedMeUp = right.fSpeedMeUp;
161}
162
164{
165 if(this == &rhs) return *this; // handle self assignment
166 //assignment operator
167 return *this;
168}
169
172{
174 fTimeStepReachedLimit = false;
175 fComputeLastPosition = false;
176 fRandomNumber = -1;
177}
178
180{
181 fpState.reset(new G4ITBrownianState());
182// G4cout << "G4DNABrownianTransportation::StartTracking : "
183// "Initialised track State" << G4endl;
186}
187
189{
190 if(verboseLevel > 0)
191 {
192 G4cout << G4endl<< GetProcessName() << ": for "
193 << setw(24) << particle.GetParticleName()
194 << "\tSubType= " << GetProcessSubType() << G4endl;
195 }
196 // Initialize water density pointer
198 GetDensityTableFor(G4Material::GetMaterial("G4_WATER"));
199
202}
203
205 const G4Step& step,
206 const G4double timeStep,
207 G4double& spaceStep)
208{
209 // G4cout << "G4ITBrownianTransportation::ComputeStep" << G4endl;
210
211 /* If this method is called, this step
212 * cannot be geometry limited.
213 * In case the step is limited by the geometry,
214 * this method should not be called.
215 * The fTransportEndPosition calculated in
216 * the method AlongStepIL should be taken
217 * into account.
218 * In order to do so, the flag IsLeadingStep
219 * is on. Meaning : this track has the minimum
220 * interaction length over all others.
221 */
222 if(GetIT(track)->GetTrackingInfo()->IsLeadingStep())
223 {
224 const G4VITProcess* ITProc = ((const G4VITProcess*) step.GetPostStepPoint()
226 G4bool makeException = true;
227
228 if(ITProc && ITProc->ProposesTimeStep()) makeException = false;
229
230 if(makeException)
231 {
232
233 G4ExceptionDescription exceptionDescription;
234 exceptionDescription << "ComputeStep is called while the track has"
235 "the minimum interaction time";
236 exceptionDescription << " so it should not recompute a timeStep ";
237 G4Exception("G4DNABrownianTransportation::ComputeStep",
238 "G4DNABrownianTransportation001", FatalErrorInArgument,
239 exceptionDescription);
240 }
241 }
242
243 State(fGeometryLimitedStep) = false;
244
245 G4Molecule* molecule = GetMolecule(track);
246
247 if(timeStep > 0)
248 {
249 spaceStep = DBL_MAX;
250
251 // TODO : generalize this process to all kind of Brownian objects
252 G4double diffCoeff = molecule->GetDiffusionCoefficient(track.GetMaterial(),
253 track.GetMaterial()->GetTemperature());
254
255 static G4double sqrt_2 = sqrt(2.);
256 G4double sqrt_Dt = sqrt(diffCoeff*timeStep);
257 G4double sqrt_2Dt = sqrt_2*sqrt_Dt;
258
259 if(State(fTimeStepReachedLimit)== true)
260 {
261 //========================================================================
262 State(fGeometryLimitedStep) = true;// important
263 //========================================================================
264 spaceStep = State(fEndPointDistance);
265 // G4cout << "State(fTimeStepReachedLimit)== true" << G4endl;
266 }
267 else
268 {
269 G4double x = G4RandGauss::shoot(0,sqrt_2Dt);
270 G4double y = G4RandGauss::shoot(0,sqrt_2Dt);
271 G4double z = G4RandGauss::shoot(0,sqrt_2Dt);
272
273 spaceStep = sqrt(x*x + y*y + z*z);
274
275 if(spaceStep >= State(fEndPointDistance))
276 {
277 //G4cout << "spaceStep >= State(fEndPointDistance)" << G4endl;
278 //======================================================================
279 State(fGeometryLimitedStep) = true;// important
280 //======================================================================
281/*
282 if(fSpeedMeUp)
283 {
284 G4cout << "SpeedMeUp" << G4endl;
285 }
286 else
287*/
288 if(fUseSchedulerMinTimeSteps == false)// jump over barrier NOT used
289 {
290#ifdef G4VERBOSE
291 if (fVerboseLevel > 1)
292 {
294 << "G4ITBrownianTransportation::ComputeStep() : "
295 << "Step was limited to boundary"
296 << RESET_COLOR
297 << G4endl;
298 }
299#endif
300 //TODO
301 if(State(fRandomNumber)>=0) // CDF is used
302 {
303 /*
304 //=================================================================
305 State(fGeometryLimitedStep) = true;// important
306 //=================================================================
307 spaceStep = State(fEndPointDistance);
308 */
309
310 //==================================================================
311 // BE AWARE THAT THE TECHNIQUE USED BELOW IS A 1D APPROXIMATION
312 // Cumulative density function for the 3D case is not yet
313 // implemented
314 //==================================================================
315// G4cout << GREEN_ON_BLUE
316// << "G4ITBrownianTransportation::ComputeStep() : "
317// << "A random number was selected"
318// << RESET_COLOR
319// << G4endl;
320 G4double value = State(fRandomNumber)+(1-State(fRandomNumber))*G4UniformRand();
321 G4double invErfc = InvErfc(value);
322 spaceStep = invErfc*2*sqrt_Dt;
323
324 if(State(fTimeStepReachedLimit)== false)
325 {
326 //================================================================
327 State(fGeometryLimitedStep) = false;// important
328 //================================================================
329 }
330 //==================================================================
331 // DEBUG
332// if(spaceStep > State(fEndPointDistance))
333// {
334// G4cout << "value = " << value << G4endl;
335// G4cout << "invErfc = " << invErfc << G4endl;
336// G4cout << "spaceStep = " << G4BestUnit(spaceStep, "Length")
337// << G4endl;
338// G4cout << "end point distance= " << G4BestUnit(State(fEndPointDistance), "Length")
339// << G4endl;
340// }
341//
342// assert(spaceStep <= State(fEndPointDistance));
343 //==================================================================
344
345 }
346 else if(fUseMaximumTimeBeforeReachingBoundary == false) // CDF is used
347 {
348 G4double min_randomNumber = Erfc(State(fEndPointDistance)/2*sqrt_Dt);
349 G4double value = min_randomNumber+(1-min_randomNumber)*G4UniformRand();
350 G4double invErfc = InvErfc(value);
351 spaceStep = invErfc*2*sqrt_Dt;
352 if(spaceStep >= State(fEndPointDistance))
353 {
354 //================================================================
355 State(fGeometryLimitedStep) = true;// important
356 //================================================================
357 }
358 else if(State(fTimeStepReachedLimit)== false)
359 {
360 //================================================================
361 State(fGeometryLimitedStep) = false;// important
362 //================================================================
363 }
364 }
365 else // CDF is NOT used
366 {
367 //==================================================================
368 State(fGeometryLimitedStep) = true;// important
369 //==================================================================
370 spaceStep = State(fEndPointDistance);
371 //TODO
372
373 /*
374 //==================================================================
375 // 1D approximation to place the brownian between its starting point
376 // and the geometry boundary
377 //==================================================================
378 double min_randomNumber = Erfc(State(fEndPointDistance)/2*sqrt_Dt);
379 double value = State(fRandomNumber)+(1-State(fRandomNumber))*G4UniformRand();
380 double invErfc = InvErfc(value*G4UniformRand());
381 spaceStep = invErfc*2*sqrt_Dt;
382 State(fGeometryLimitedStep) = false;
383 */
384 }
385 }
386
387 State(fTransportEndPosition)= spaceStep*
388// step.GetPostStepPoint()->GetMomentumDirection()
390 + track.GetPosition();
391 }
392 else
393 {
394 //======================================================================
395 State(fGeometryLimitedStep) = false;// important
396 //======================================================================
397 State(fTransportEndPosition)= spaceStep*step.GetPostStepPoint()->
398 GetMomentumDirection() + track.GetPosition();
399 }
400 }
401 }
402 else
403 {
404 spaceStep = 0.;
405 State(fTransportEndPosition) = track.GetPosition();
406 State(fGeometryLimitedStep) = false;
407 }
408
409 State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime()
410 + timeStep;
411 State(fEndGlobalTimeComputed) = true;
412
413#ifdef G4VERBOSE
414 // DEBUG
415 if (fVerboseLevel > 1)
416 {
418 << "G4ITBrownianTransportation::ComputeStep() : "
419 << " trackID : " << track.GetTrackID() << " : Molecule name: "
420 << molecule->GetName() << G4endl
421 << "Initial position:" << G4BestUnit(track.GetPosition(), "Length")
422 << G4endl
423 << "Initial direction:" << track.GetMomentumDirection() << G4endl
424 << "Final position:" << G4BestUnit(State(fTransportEndPosition), "Length")
425 << G4endl
426 << "Initial magnitude:" << G4BestUnit(track.GetPosition().mag(), "Length")
427 << G4endl
428 << "Final magnitude:" << G4BestUnit(State(fTransportEndPosition).mag(), "Length")
429 << G4endl
430 << "Diffusion length : "
431 << G4BestUnit(spaceStep, "Length")
432 << " within time step : " << G4BestUnit(timeStep,"Time")
433 << G4endl
434 << "State(fTimeStepReachedLimit)= " << State(fTimeStepReachedLimit) << G4endl
435 << "State(fGeometryLimitedStep)=" << State(fGeometryLimitedStep) << G4endl
436 << "End point distance was: " << G4BestUnit(State(fEndPointDistance), "Length")
437 << G4endl
438 << RESET_COLOR
439 << G4endl<< G4endl;
440 }
441#endif
442
443//==============================================================================
444// DEBUG
445//assert(spaceStep < State(fEndPointDistance)
446// || (spaceStep >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
447//assert(track.GetMomentumDirection() == State(fTransportEndMomentumDir));
448//==============================================================================
449}
450
452 const G4Step& step)
453{
455
456#ifdef G4VERBOSE
457 // DEBUG
458 if (fVerboseLevel > 1)
459 {
460 G4cout << GREEN_ON_BLUE << "G4ITBrownianTransportation::PostStepDoIt() :"
461 << " trackID : " << track.GetTrackID() << " Molecule name: "
462 << GetMolecule(track)->GetName() << G4endl;
463 G4cout << "Diffusion length : "
464 << G4BestUnit(step.GetStepLength(), "Length")
465 <<" within time step : " << G4BestUnit(step.GetDeltaTime(),"Time")
466 << "\t Current global time : "
467 << G4BestUnit(track.GetGlobalTime(),"Time")
468 << RESET_COLOR
469 << G4endl<< G4endl;
470 }
471#endif
472 return &fParticleChange;
473}
474
476{
477
478#ifdef DEBUG_MEM
479 MemStat mem_first = MemoryUsage();
480#endif
481
482#ifdef G4VERBOSE
483 // DEBUG
484 if (fVerboseLevel > 1)
485 {
486 G4cout << GREEN_ON_BLUE << setw(18)
487 << "G4DNABrownianTransportation::Diffusion :" << setw(8)
488 << GetIT(track)->GetName() << "\t trackID:" << track.GetTrackID()
489 << "\t" << " Global Time = "
490 << G4BestUnit(track.GetGlobalTime(), "Time")
491 << RESET_COLOR
492 << G4endl
493 << G4endl;
494 }
495#endif
496
497/*
498 fParticleChange.ProposePosition(State(fTransportEndPosition));
499 //fParticleChange.ProposeEnergy(State(fTransportEndKineticEnergy));
500 fParticleChange.SetMomentumChanged(State(fMomentumChanged));
501
502 fParticleChange.ProposeGlobalTime(State(fCandidateEndGlobalTime));
503 fParticleChange.ProposeLocalTime(State(fCandidateEndGlobalTime));
504 fParticleChange.ProposeTrueStepLength(track.GetStepLength());
505*/
507
508 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
509
510 if (waterDensity == 0.0)
511 {
513 {
514 // Let the user Brownian action class decide what to do
517 return;
518 }
519 else
520 {
521#ifdef G4VERBOSE
522 if(fVerboseLevel)
523 {
524 G4cout << "A track is outside water material : trackID = "
525 << track.GetTrackID() << " (" << GetMolecule(track)->GetName() <<")"
526 << G4endl;
527 G4cout << "Local Time : " << G4BestUnit(track.GetGlobalTime(), "Time")
528 << G4endl;
529 G4cout << "Step Number :" << track.GetCurrentStepNumber() << G4endl;
530 }
531#endif
534 return;// &fParticleChange is the final returned object
535 }
536 }
537
538
539 #ifdef DEBUG_MEM
540 MemStat mem_intermediaire = MemoryUsage();
541 mem_diff = mem_intermediaire-mem_first;
542 G4cout << "\t\t\t >> || MEM || In G4DNABrownianTransportation::Diffusion "
543 "after dealing with waterDensity for "<< track.GetTrackID()
544 << ", diff is : " << mem_diff << G4endl;
545 #endif
546
548 State(fMomentumChanged) = true;
550 //
551 #ifdef DEBUG_MEM
552 mem_intermediaire = MemoryUsage();
553 mem_diff = mem_intermediaire-mem_first;
554 G4cout << "\t\t\t >> || MEM || In G4DNABrownianTransportation::"
555 "After proposing new direction to fParticleChange for "
556 << track.GetTrackID() << ", diff is : " << mem_diff << G4endl;
557 #endif
558
559 return;// &fParticleChange is the final returned object
560}
561
562// NOT USED
564 G4double& presafety,
565 G4double limit)
566{
567 G4double res = DBL_MAX;
568 if(track.GetVolume() != fpSafetyHelper->GetWorldVolume())
569 {
570 G4TrackStateManager& trackStateMan = GetIT(track)->GetTrackingInfo()
572 fpSafetyHelper->LoadTrackState(trackStateMan);
574 track.GetStep()->GetPreStepPoint()->GetPosition(),
575 track.GetMomentumDirection(),
576 limit, presafety);
578 }
579 return res;
580}
581
583 G4double previousStepSize,
584 G4double currentMinimumStep,
585 G4double& currentSafety,
586 G4GPILSelection* selection)
587{
588#ifdef G4VERBOSE
589 if(fVerboseLevel)
590 {
591 G4cout << " G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength - track ID: "
592 << track.GetTrackID() << G4endl;
593 G4cout << "In volume : " << track.GetVolume()->GetName()
594 << " position : " << G4BestUnit(track.GetPosition() , "Length") << G4endl;
595 }
596#endif
597
598 G4double geometryStepLength =
600 track, previousStepSize, currentMinimumStep, currentSafety,
601 selection);
602
603 if(geometryStepLength==0)
604 {
605// G4cout << "geometryStepLength==0" << G4endl;
606 if(State(fGeometryLimitedStep))
607 {
608// G4cout << "if(State(fGeometryLimitedStep))" << G4endl;
609 G4TouchableHandle newTouchable = new G4TouchableHistory;
610
611 newTouchable->UpdateYourself(State(fCurrentTouchableHandle)->GetVolume(),
612 State(fCurrentTouchableHandle)->GetHistory());
613
614 fLinearNavigator->SetGeometricallyLimitedStep();
615 fLinearNavigator->LocateGlobalPointAndUpdateTouchableHandle(
616 track.GetPosition(), track.GetMomentumDirection(),
617 newTouchable, true);
618
619 if(newTouchable->GetVolume() == 0)
620 {
621 return 0;
622 }
623
624 State(fCurrentTouchableHandle) = newTouchable;
625
626 //=======================================================================
627 // TODO: speed up navigation update
628// geometryStepLength = fLinearNavigator->ComputeStep(track.GetPosition(),
629// track.GetMomentumDirection(),
630// currentMinimumStep,
631// currentSafety);
632 //=======================================================================
633
634
635 //=======================================
636 // Longer but safer ...
637 geometryStepLength =
639 track, previousStepSize, currentMinimumStep, currentSafety,
640 selection);
641
642 }
643 }
644
645 //============================================================================
646 // DEBUG
647 // G4cout << "geometryStepLength: " << G4BestUnit(geometryStepLength, "Length")
648 // << " | trackID: " << track.GetTrackID()
649 // << G4endl;
650 //============================================================================
651
652 G4double diffusionCoefficient = 0;
653
654 /*
655 G4Material* material = track.GetMaterial();
656 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
657
658 if (waterDensity == 0.0)
659 {
660 if(fpBrownianAction)
661 {
662 diffusionCoefficient = fpBrownianAction->GetDiffusionCoefficient(material,
663 GetMolecule(track));
664 }
665
666 if(diffusionCoefficient <= 0)
667 {
668 State(fGeometryLimitedStep) = false;
669 State(theInteractionTimeLeft) = 0;
670 State(fTransportEndPosition) = track.GetPosition();
671 return 0;
672 }
673
674 }
675 else
676 {
677 diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
678 }
679 */
680
681 diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
682
683 // To avoid divide by zero of diffusionCoefficient
684 if(diffusionCoefficient <= 0)
685 {
686 State(fGeometryLimitedStep) = false;
688 State(fTransportEndPosition) = track.GetPosition();
689 return 0;
690 }
691
692
693 State(fComputeLastPosition) = false;
694 State(fTimeStepReachedLimit) = false;
695
696 if (State(fGeometryLimitedStep))
697 {
698 // 95 % of the space step distribution is lower than
699 // d_95 = 2 * sqrt(2*D*t)
700 // where t is the corresponding time step
701 // so by inversion :
703 {
704 if(fSpeedMeUp)
705 {
706 State(theInteractionTimeLeft) = (geometryStepLength * geometryStepLength)
707 / (diffusionCoefficient); // d_50 - use straight line
708 }
709 else
710 {
711 State(theInteractionTimeLeft) = (currentSafety * currentSafety)
712 / (diffusionCoefficient); // d_50 - use safety
713
714 //=====================================================================
715 // State(theInteractionTimeLeft) = (currentSafety * currentSafety)
716 // / (8 * diffusionCoefficient); // d_95
717 //=====================================================================
718 }
719 State(fComputeLastPosition) = true;
720 }
721 else
722 // Will use a random time - this is precise but long to compute in certain
723 // circumstances (many particles - small volumes)
724 {
725 State(fRandomNumber) = G4UniformRand();
726 State(theInteractionTimeLeft) = 1 / (4 * diffusionCoefficient)
727 * pow(geometryStepLength / InvErfc(State(fRandomNumber)),2);
728
729 State(fTransportEndPosition) = geometryStepLength*
730 track.GetMomentumDirection() + track.GetPosition();
731 }
732
734 {
735 G4double minTimeStepAllowed = G4VScheduler::Instance()->GetLimitingTimeStep();
736 //========================================================================
737 // TODO
738// double currentMinTimeStep = G4VScheduler::Instance()->GetTimeStep();
739 //========================================================================
740
741 if (State(theInteractionTimeLeft) < minTimeStepAllowed)
742 {
743 State(theInteractionTimeLeft) = minTimeStepAllowed;
744 State(fTimeStepReachedLimit) = true;
745 State(fComputeLastPosition) = true;
746 }
747 }
749 // TODO: find a better way when fForceLimitOnMinTimeSteps is not used
750 {
751 State(fTimeStepReachedLimit) = true;
754 {
755 State(fComputeLastPosition) = true;
756 }
757 }
758
759 State(fCandidateEndGlobalTime) =
761
762 State(fEndGlobalTimeComputed) = true; // MK: ADDED ON 20/11/2014
763
764 State(fPathLengthWasCorrected) = false;
765 }
766 else
767 {
768 // Transform geometrical step
769 geometryStepLength = 2
770 * sqrt(diffusionCoefficient * State(theInteractionTimeLeft))
772 State(fPathLengthWasCorrected) = true;
773 //State(fEndPointDistance) = geometryStepLength;
774 State(fTransportEndPosition) = geometryStepLength*
775 track.GetMomentumDirection() + track.GetPosition();
776 }
777
778#ifdef G4VERBOSE
779 // DEBUG
780 if (fVerboseLevel > 1)
781 {
783 << "G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength = "
784 << G4BestUnit(geometryStepLength, "Length")
785 << " | trackID = "
786 << track.GetTrackID()
787 << RESET_COLOR
788 << G4endl;
789 }
790#endif
791
792// assert(geometryStepLength < State(fEndPointDistance)
793// || (geometryStepLength >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
794
795 return geometryStepLength;
796}
797
799//
800// Initialize ParticleChange (by setting all its members equal
801// to corresponding members in G4Track)
804 const G4Step& step)
805{
806#ifdef DEBUG_MEM
807 MemStat mem_first, mem_second, mem_diff;
808#endif
809
810#ifdef DEBUG_MEM
811 mem_first = MemoryUsage();
812#endif
813
814 if (GetIT(track)->GetTrackingInfo()->IsLeadingStep()
815 && State(fComputeLastPosition)
816 && State(fGeometryLimitedStep)//Hoang added 26/8/2019
817 )
818 {
819 //==========================================================================
820 // DEBUG
821 //
822// assert(fabs(State(theInteractionTimeLeft)-
823// G4VScheduler::Instance()->GetTimeStep()) < DBL_EPSILON);
824 //==========================================================================
825
826 G4double spaceStep = DBL_MAX;
827
829 {
830 spaceStep = State(fEndPointDistance);
831 State(fGeometryLimitedStep) = true;
832 }
833 else
834 {
835 G4double diffusionCoefficient = GetMolecule(track)->
836 GetDiffusionCoefficient();
837
838 G4double sqrt_2Dt= sqrt(2 * diffusionCoefficient * State(theInteractionTimeLeft));
839 G4double x = G4RandGauss::shoot(0, sqrt_2Dt);
840 G4double y = G4RandGauss::shoot(0, sqrt_2Dt);
841 G4double z = G4RandGauss::shoot(0, sqrt_2Dt);
842
843 spaceStep = sqrt(x * x + y * y + z * z);
844
845 if(spaceStep >= State(fEndPointDistance))
846 {
847 State(fGeometryLimitedStep) = true;
848 if (
849 //fSpeedMeUp == false&&
851 && spaceStep >= State(fEndPointDistance))
852 {
853 spaceStep = State(fEndPointDistance);
854 }
855 }
856 else
857 {
858 State(fGeometryLimitedStep) = false;
859 }
860 }
861
862// assert( (spaceStep < State(fEndPointDistance) && State(fGeometryLimitedStep) == false)
863//|| (spaceStep >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
864
865 // Calculate final position
866 //
867 State(fTransportEndPosition) = track.GetPosition()
868 + spaceStep * track.GetMomentumDirection();
869 }
870
871 if(fVerboseLevel)
872 {
874 << "G4DNABrownianTransportation::AlongStepDoIt: "
875 "GeometryLimitedStep = "
876 << State(fGeometryLimitedStep)
877 << RESET_COLOR
878 << G4endl;
879 }
880
881// static G4ThreadLocal G4int noCalls = 0;
882// noCalls++;
883
884// fParticleChange.Initialize(track);
885
887
888#ifdef DEBUG_MEM
889 MemStat mem_intermediaire = MemoryUsage();
890 mem_diff = mem_intermediaire-mem_first;
891 G4cout << "\t\t\t >> || MEM || After calling G4ITTransportation::"
892 "AlongStepDoIt for "<< track.GetTrackID() << ", diff is : "
893 << mem_diff << G4endl;
894#endif
895
896 if(track.GetStepLength() != 0 // && State(fGeometryLimitedStep)
897 //========================================================================
898 // TODO: avoid changing direction after too small time steps
899// && (G4VScheduler::Instance()->GetTimeStep() > fInternalMinTimeStep
900// || fSpeedMeUp == false)
901 //========================================================================
902 )
903 {
904 Diffusion(track);
905 }
906 //else
907 //{
908 // fParticleChange.ProposeMomentumDirection(State(fTransportEndMomentumDir));
909 //}
910/*
911 if (State(fParticleIsLooping))
912 {
913 if ((State(fNoLooperTrials)>= fThresholdTrials))
914 {
915 fParticleChange.ProposeTrackStatus(fStopAndKill);
916 State(fNoLooperTrials) = 0;
917#ifdef G4VERBOSE
918 if ((fVerboseLevel > 1))
919 {
920 G4cout
921 << " G4DNABrownianTransportation is killing track that is looping or stuck "
922 << G4endl;
923 G4cout << " Number of trials = " << State(fNoLooperTrials)
924 << " No of calls to AlongStepDoIt = " << noCalls
925 << G4endl;
926 }
927#endif
928 }
929 else
930 {
931 State(fNoLooperTrials)++;
932 }
933 }
934 else
935 {
936 State(fNoLooperTrials)=0;
937 }
938*/
939#ifdef DEBUG_MEM
940 mem_intermediaire = MemoryUsage();
941 mem_diff = mem_intermediaire-mem_first;
942 G4cout << "\t\t\t >> || MEM || After calling G4DNABrownianTransportation::"
943 "Diffusion for "<< track.GetTrackID() << ", diff is : "
944 << mem_diff << G4endl;
945#endif
946
947 return &fParticleChange;
948}
949
static G4double InvErf(G4double x)
#define GREEN_ON_BLUE
static G4double Erfc(G4double x)
#define State(theXInfo)
static G4double InvErfc(G4double x)
#define RESET_COLOR
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4GPILSelection
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
@ fLowEnergyBrownianTransportation
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
G4ThreeVector G4RandomDirection()
static constexpr double picosecond
Definition: G4SIunits.hh:141
#define G4BestUnit(a, b)
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double mag() const
static double erf(double x)
static double inverseErf(double t)
virtual void Transport(const G4Track &, G4ParticleChangeForTransport &)=0
virtual void ComputeStep(const G4Track &, const G4Step &, const G4double, G4double &)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &)
G4DNABrownianTransportation & operator=(const G4DNABrownianTransportation &)
void Diffusion(const G4Track &track)
G4DNABrownianTransportation(const G4String &aName="DNABrownianTransportation", G4int verbosityLevel=0)
const std::vector< G4double > * fpWaterDensity
G4double ComputeGeomLimit(const G4Track &track, G4double &presafety, G4double limit)
virtual void StartTracking(G4Track *aTrack)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
static G4DNAMolecularMaterial * Instance()
G4VPhysicalVolume * GetWorldVolume()
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4ParticleChangeForTransport fParticleChange
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetInstantiateProcessState(G4bool flag)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4ITSafetyHelper * fpSafetyHelper
G4ITNavigator * fLinearNavigator
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
virtual void StartTracking(G4Track *aTrack)
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:143
virtual const G4String & GetName() const =0
G4double GetTemperature() const
Definition: G4Material.hh:178
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:686
const G4String & GetName() const
Definition: G4Molecule.cc:338
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:516
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
void SetMomentumChanged(G4bool b)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:62
G4double GetDeltaTime() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
virtual void LoadTrackState(G4TrackStateManager &manager)
virtual void ResetTrackState()
G4int GetTrackID() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4int GetCurrentStepNumber() const
G4Material * GetMaterial() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetStepLength() const
const G4Step * GetStep() const
G4TrackStateManager & GetTrackStateManager()
G4bool ProposesTimeStep() const
G4double * theInteractionTimeLeft
G4shared_ptr< G4ProcessState > fpState
void ProposeTrackStatus(G4TrackStatus status)
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
static G4VScheduler * Instance()
Definition: G4VScheduler.cc:48
virtual double GetLimitingTimeStep() const
Definition: G4VScheduler.hh:84
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
virtual void UpdateYourself(G4VPhysicalVolume *pPhysVol, const G4NavigationHistory *history=nullptr)
Definition: G4VTouchable.cc:66
ThreeVector shoot(const G4int Ap, const G4int Af)
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62