G4DNABrownianTransportation.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4DNABrownianTransportation.cc 64374 2012-10-31 16:37:23Z gcosmo $
00027 //
00028 // Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr) 
00029 //
00030 // WARNING : This class is released as a prototype.
00031 // It might strongly evolve or even disapear in the next releases.
00032 //
00033 // History:
00034 // -----------
00035 // 10 Oct 2011 M.Karamitros created
00036 //
00037 // -------------------------------------------------------------------
00038 
00041 
00042 #include <CLHEP/Random/Stat.h>
00043 
00044 #include "G4DNABrownianTransportation.hh"
00045 #include "G4PhysicalConstants.hh"
00046 #include "G4SystemOfUnits.hh"
00047 #include "Randomize.hh"
00048 #include "G4Molecule.hh"
00049 #include "G4RandomDirection.hh"
00050 #include "G4ParticleTable.hh"
00051 #include "G4SafetyHelper.hh"
00052 #include "G4TransportationManager.hh"
00053 #include "G4UnitsTable.hh"
00054 #include "G4NistManager.hh"
00055 #include "G4DNAMolecularMaterial.hh"
00056 
00057 using namespace std;
00058 
00059 #ifndef State
00060 #define State(theXInfo) (fpBrownianState->theXInfo)
00061 #endif
00062 
00063 //#ifndef State
00064 //#define State(theXInfo) (fTransportationState->theXInfo)
00065 //#endif
00066 
00067 //COLOR FOR DEBUGING
00068 //#define RED_ON_WHITE  "\033[0;31m"
00069 //#define GREEN "\033[32;40m"
00070 #define GREEN_ON_BLUE "\033[1;32;44m"
00071 #define RESET "\033[0m"
00072 
00073 static double InvErf(double x)
00074 {
00075     return CLHEP::HepStat::inverseErf(x);
00076 }
00077 
00078 static double InvErfc(double x)
00079 {
00080     return CLHEP::HepStat::inverseErf(1.-x);
00081 }
00082 
00083 G4DNABrownianTransportation::G4DNABrownianTransportation(const G4String& aName, G4int verbosity) :
00084     G4ITTransportation(aName, verbosity),
00085     InitProcessState(fpBrownianState, fTransportationState)
00086 {
00087     //ctor
00088     SetProcessSubType(61);
00089     verboseLevel = 1;
00090     fUseMaximumTimeBeforeReachingBoundary = true;
00091     fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
00092     fpWaterDensity  = 0;
00093 }
00094 
00095 G4DNABrownianTransportation::~G4DNABrownianTransportation()
00096 {;}
00097 
00098 G4DNABrownianTransportation::G4DNABrownianTransportation(const G4DNABrownianTransportation& right) :
00099     G4ITTransportation(right),
00100     InitProcessState(fpBrownianState, fTransportationState)
00101 {
00102     //copy ctor
00103     SetProcessSubType(61);
00104     verboseLevel = right.verboseLevel;
00105     fUseMaximumTimeBeforeReachingBoundary = right.fUseMaximumTimeBeforeReachingBoundary;
00106     fNistWater = right.fNistWater;
00107     fpWaterDensity  = right.fpWaterDensity ;
00108 }
00109 
00110 G4DNABrownianTransportation& G4DNABrownianTransportation::operator=(const G4DNABrownianTransportation& rhs)
00111 {
00112     if (this == &rhs) return *this; // handle self assignment
00113     //assignment operator
00114     return *this;
00115 }
00116 
00117 G4DNABrownianTransportation::G4ITBrownianState::G4ITBrownianState() : G4ITTransportationState()
00118 {
00119     fPathLengthWasCorrected = false;
00120 }
00121 
00122 void G4DNABrownianTransportation::StartTracking(G4Track* track)
00123 {
00124     fpState = new G4ITBrownianState();
00125     SetInstantiateProcessState(false);
00126     G4ITTransportation::StartTracking(track);
00127 }
00128 
00129 void G4DNABrownianTransportation::BuildPhysicsTable(const G4ParticleDefinition& particle)
00130 {
00131     if(verboseLevel > 0)
00132     {
00133         G4cout << G4endl << GetProcessName() << ":   for  "
00134                << setw(24) << particle.GetParticleName()
00135                << "\tSubType= " << GetProcessSubType()   << G4endl;
00136     }
00137 
00138     // Initialize water density pointer
00139     fpWaterDensity = G4DNAMolecularMaterial::Instance()->GetDensityTableFor(G4Material::GetMaterial("G4_WATER"));
00140 }
00141 
00142 void G4DNABrownianTransportation::ComputeStep(const G4Track& track,
00143                                               const G4Step& step,
00144                                               const double timeStep,
00145                                               double& spaceStep)
00146 {
00147     // G4cout << "G4ITBrownianTransportation::ComputeStep" << G4endl;
00148 
00149     // If this method is called, this step
00150     // cannot be geometry limited.
00151     // In case the step is limited by the geometry,
00152     // this method should not be called.
00153     // The fTransportEndPosition calculated in
00154     // the method AlongStepIL should be taken
00155     // into account.
00156     // In order to do so, the flag IsLeadingStep
00157     // is on. Meaning : this track has the minimum
00158     // interaction length over all others.
00159     if(GetIT(track)->GetTrackingInfo()->IsLeadingStep())
00160     {
00161         const G4VITProcess* ITProc = ((const G4VITProcess*) step.GetPostStepPoint()->GetProcessDefinedStep());
00162         bool makeException = true;
00163 
00164         if(ITProc && ITProc->ProposesTimeStep()) makeException = false;
00165 
00166         if(makeException)
00167         {
00168 
00169         G4ExceptionDescription exceptionDescription ;
00170         exceptionDescription << "ComputeStep is called while the track has the minimum interaction time";
00171         exceptionDescription << " so it should not recompute a timeStep ";
00172         G4Exception("G4DNABrownianTransportation::ComputeStep","G4DNABrownianTransportation001",
00173                     FatalErrorInArgument,exceptionDescription);
00174         }
00175     }
00176 
00177     State(fGeometryLimitedStep) = false;
00178     // TODO : generalize this process to all kind of brownian objects
00179     //    G4ITBrownianObject* ITBrown = GetITBrownianObject(track) ;
00180     //    G4double diffCoeff = ITBrown->GetDiffusionCoefficient(track.GetMaterial());
00181     G4Molecule* molecule = GetMolecule(track) ;
00182     G4double diffCoeff = molecule->GetDiffusionCoefficient();
00183 
00184     if(timeStep > 0)
00185     {
00186         spaceStep = DBL_MAX;
00187 
00188         while(spaceStep > State(endpointDistance))
00189             // Probably inefficient when the track is close close to boundaries
00190             // it goes with fUserMaximumTimeBeforeReachingBoundary == false
00191             // fUserMaximumTimeBeforeReachingBoundary == true, it should never loop
00192         {
00193             G4double x = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
00194             G4double y = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
00195             G4double z = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
00196 
00197             spaceStep = sqrt(x*x + y*y + z*z);
00198         }
00199         //        State(fTransportEndPosition).set(x + track.GetPosition().x(),
00200         //                                         y + track.GetPosition().y(),
00201         //                                         z + track.GetPosition().z());
00202 
00203         State(fTransportEndPosition)= spaceStep*step.GetPostStepPoint()->GetMomentumDirection() + track.GetPosition();
00204     }
00205     else
00206     {
00207         spaceStep = 0. ;
00208         State(fTransportEndPosition) =  track.GetPosition() ;
00209     }
00210 
00211     State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime() + timeStep ;
00212     State(fEndGlobalTimeComputed) = true ;
00213 
00214 #ifdef G4VERBOSE
00215     //    DEBUG
00216     if(fVerboseLevel>1)
00217     {
00218         G4cout<< GREEN_ON_BLUE
00219               << "G4ITBrownianTransportation::ComputeStep() : "
00220               << " trackID : "         << track.GetTrackID()
00221               << " : Molecule name: "  << molecule-> GetName()
00222               << G4endl
00223               << "Diffusion length : " << G4BestUnit(spaceStep, "Length")
00224               << " within time step : "  << G4BestUnit(timeStep,"Time")
00225               << RESET
00226               << G4endl<< G4endl;
00227     }
00228 #endif
00229 }
00230 
00231 G4VParticleChange* G4DNABrownianTransportation::PostStepDoIt( const G4Track& track, const G4Step& step)
00232 {
00233     G4ITTransportation::PostStepDoIt(track,step);
00234 
00235 #ifdef G4VERBOSE
00236     //    DEBUG
00237     if(fVerboseLevel>1)
00238     {
00239         G4cout<< GREEN_ON_BLUE
00240               << "G4ITBrownianTransportation::PostStepDoIt() :"
00241               << " trackID : "        << track.GetTrackID()
00242               << " Molecule name: "   << GetMolecule(track)-> GetName() << G4endl;
00243         G4cout<< "Diffusion length : "<< G4BestUnit(step.GetStepLength(),"Length") <<" within time step : "
00244               << G4BestUnit(step.GetDeltaTime(),"Time") << "\t"
00245               << " Current global time : " << G4BestUnit(track.GetGlobalTime(),"Time")
00246               << RESET
00247               << G4endl<< G4endl;
00248     }
00249 #endif
00250     return &fParticleChange ;
00251 }
00252 
00253 void G4DNABrownianTransportation::Diffusion(
00254         const G4Track& track)
00255 {
00256 
00257 #ifdef G4VERBOSE
00258     // DEBUG
00259     if (fVerboseLevel>1)
00260     {
00261         G4cout<< GREEN_ON_BLUE
00262               << setw(18)<< "G4DNABrownianTransportation::Diffusion :"
00263               << setw(8) <<  GetIT(track)->GetName()
00264               << "\t trackID:"      << track.GetTrackID() <<"\t"
00265               << " Global Time = "  << G4BestUnit(track.GetGlobalTime(),"Time")
00266               << RESET
00267               << G4endl<< G4endl;
00268     }
00269 #endif
00270 
00271     G4Material* material = track.GetMaterial();
00272 //    if (material != fNistWater && material->GetBaseMaterial() != fNistWater)
00273 
00274    G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
00275 
00276    if(waterDensity == 0.0)
00277 //  if (material == nistwater || material->GetBaseMaterial() == nistwater)
00278     {
00279         G4cout << "A track is outside water material : trackID"<< track.GetTrackID() << " (" << GetMolecule(track)->GetName()  <<")" << G4endl;
00280         G4cout << "Local Time : "<<  (track.GetLocalTime()) /s<<G4endl;
00281         G4cout<< "Step Number :" << track.GetCurrentStepNumber() <<G4endl;
00282 
00283         fParticleChange.ProposeEnergy(0.) ;
00284         fParticleChange.ProposeTrackStatus(fStopButAlive);
00285 
00286         // Got pb with :
00287         // fParticleChange.ProposeTrackStatus(fStopAndKill);
00288         // It makes the tracks going straight without killing them
00289 
00290         return ; // &fParticleChange is the final returned object
00291     }
00292 
00293     G4double costheta = (2*G4UniformRand()-1);
00294     G4double theta = acos (costheta);
00295     G4double phi = 2*pi*G4UniformRand();
00296 
00297     G4double xMomentum = cos(phi)* sin(theta);
00298     G4double yMomentum = sin(theta)*sin(phi);
00299     G4double zMomentum = costheta;
00300 
00301     fParticleChange.ProposeMomentumDirection(xMomentum, yMomentum, zMomentum);
00302     State(fMomentumChanged) = true;
00303     fParticleChange.SetMomentumChanged(true) ;
00304 
00305     //    G4cout << "BROWN : Propose new direction :" << G4ThreeVector(xMomentum, yMomentum, zMomentum) << G4endl;
00306 
00307     // Alternative
00308     //fParticleChange.ProposeMomentumDirection(G4RandomDirection());
00309 
00310     return; // &fParticleChange is the final returned object
00311 }
00312 
00313 
00314 G4double G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(
00315         const G4Track& track,
00316         G4double previousStepSize,
00317         G4double currentMinimumStep,
00318         G4double& currentSafety,
00319         G4GPILSelection* selection)
00320 {
00321     G4double geometryStepLength = G4ITTransportation::AlongStepGetPhysicalInteractionLength(track,
00322                                                                                             previousStepSize,
00323                                                                                             currentMinimumStep,
00324                                                                                             currentSafety,
00325                                                                                             selection);
00326 
00327     G4double diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
00328     //    G4double diffusionCoefficient = GetITBrownianObject(track)->GetDiffusionCoefficient(track.GetMaterial());
00329 
00330     if(State(fGeometryLimitedStep))
00331     {
00332         // 99 % of the space step distribution is lower than
00333         // d_99 = 8 * sqrt(D*t)
00334         // where t is the corresponding time step
00335         // so by inversion :
00336         if(fUseMaximumTimeBeforeReachingBoundary)
00337         {
00338             State(theInteractionTimeLeft)     = (geometryStepLength*geometryStepLength)/(64 * diffusionCoefficient);
00339         }
00340         else // Will use a random time
00341         {
00342             State(theInteractionTimeLeft) = 1/(4*diffusionCoefficient) * pow(geometryStepLength/InvErfc(G4UniformRand()),2);
00343         }
00344 
00345         State(fCandidateEndGlobalTime) = track.GetGlobalTime() + State(theInteractionTimeLeft);
00346         State(fPathLengthWasCorrected) = false;
00347     }
00348     else
00349     {
00350         geometryStepLength  = 2*sqrt(diffusionCoefficient*State(theInteractionTimeLeft) ) *InvErf(G4UniformRand());
00351         State(fPathLengthWasCorrected) = true;
00352         State(endpointDistance) = geometryStepLength;
00353     }
00354 
00355     return geometryStepLength ;
00356 }
00357 
00359 //
00360 //   Initialize ParticleChange  (by setting all its members equal
00361 //                               to corresponding members in G4Track)
00362 G4VParticleChange* G4DNABrownianTransportation::AlongStepDoIt( const G4Track& track,
00363                                                                const G4Step&  step )
00364 {
00365     G4ITTransportation::AlongStepDoIt(track,step);
00366     Diffusion(track);
00367     return &fParticleChange;
00368 }

Generated on Mon May 27 17:48:02 2013 for Geant4 by  doxygen 1.4.7