G4DNAMolecularReaction.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: G4DNAMolecularReaction.cc 64057 2012-10-30 15:04:49Z gcosmo $
00027 //
00028 // Author: Mathieu Karamitros (kara@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 
00039 #include "G4DNAMolecularReaction.hh"
00040 #include "G4DNAMolecularReactionTable.hh"
00041 #include "G4VDNAReactionModel.hh"
00042 #include "G4ITManager.hh"
00043 #include "G4UnitsTable.hh"
00044 
00045 using namespace std;
00046 
00047 G4DNAMolecularReaction::G4DNAMolecularReaction():G4VITReactionProcess(),
00048     fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
00049 {
00050     //ctor
00051     fVerbose        = 0;
00052     fReactionModel  = 0;
00053     fReactionRadius = -1;
00054     fDistance       = -1;
00055 }
00056 
00057 G4DNAMolecularReaction::~G4DNAMolecularReaction()
00058 {
00059     //dtor
00060     if(fpChanges) delete fpChanges;
00061 }
00062 
00063 G4DNAMolecularReaction::G4DNAMolecularReaction(const G4DNAMolecularReaction& other):G4VITReactionProcess(other),
00064     fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
00065 {
00066     //copy ctor
00067     fVerbose            = other.fVerbose ;
00068     fMolReactionTable   = other.fMolReactionTable;
00069     fReactionModel      = 0;
00070     fReactionRadius     = -1;
00071     fDistance           = -1;
00072 }
00073 
00074 G4DNAMolecularReaction& G4DNAMolecularReaction::operator=(const G4DNAMolecularReaction& rhs)
00075 {
00076     if (this == &rhs) return *this; // handle self assignment
00077 
00078     fVerbose            = rhs.fVerbose ;
00079     fMolReactionTable   = rhs.fMolReactionTable;
00080     fReactionRadius     = -1;
00081     fDistance           = -1;
00082 
00083     //assignment operator
00084     return *this;
00085 }
00086 
00087 G4bool G4DNAMolecularReaction::TestReactibility(const G4Track& trackA,
00088                                                 const G4Track& trackB,
00089                                                 const double currentStepTime,
00090                                                 const double /*previousStepTime*/,
00091                                                 bool userStepTimeLimit) /*const*/
00092 {
00093     G4Molecule* moleculeA = GetMolecule(trackA);
00094     G4Molecule* moleculeB = GetMolecule(trackB);
00095 
00096     if(!fReactionModel)
00097     {
00098         G4ExceptionDescription exceptionDescription ("You have to give a reaction model to the molecular reaction process");
00099         G4Exception("G4DNAMolecularReaction::TestReactibility","MolecularReaction001",
00100                     FatalErrorInArgument,exceptionDescription);
00101         return false; // makes coverity happy
00102     }
00103     if(fMolReactionTable==0)
00104     {
00105         G4ExceptionDescription exceptionDescription ("You have to give a reaction table to the molecular reaction process");
00106         G4Exception("G4DNAMolecularReaction::TestReactibility","MolecularReaction002",
00107                     FatalErrorInArgument,exceptionDescription);
00108         return false; // makes coverity happy
00109     }
00110 
00111     // Retrieve reaction range
00112     fReactionRadius = -1 ; // reaction Range
00113     fReactionRadius = fReactionModel -> GetReactionRadius(moleculeA, moleculeB);
00114 
00115     fDistance = -1 ; // separation distance
00116 
00117     if(currentStepTime == 0.)
00118     {
00119         userStepTimeLimit = false;
00120     }
00121 
00122     G4bool output = fReactionModel->FindReaction(trackA, trackB, fReactionRadius,fDistance, userStepTimeLimit);
00123 
00124 #ifdef G4VERBOSE
00125     //    DEBUG
00126     if(fVerbose > 1)
00127     {
00128         G4cout << "\033[1;39;36m" << "G4MolecularReaction "<< G4endl;
00129         G4cout << "FindReaction returned : " << G4BestUnit(output,"Length") << G4endl;
00130         G4cout << " reaction radius : " << G4BestUnit(fReactionRadius,"Length")
00131                << " real distance : " << G4BestUnit((trackA.GetPosition() - trackB.GetPosition()).mag(), "Length")
00132                << " calculated distance by model (= -1 if real distance > reaction radius and the user limitation step is not reached) : "
00133                << G4BestUnit(fDistance,"Length")
00134                << G4endl;
00135 
00136         G4cout << "TrackID A : " << trackA.GetTrackID()
00137                << ", TrackID B : " << trackB.GetTrackID()
00138                << " | MolA " << moleculeA->GetName()
00139                << ", MolB " << moleculeB->GetName()
00140                <<"\033[0m\n"
00141                << G4endl;
00142         G4cout << "--------------------------------------------" << G4endl;
00143     }
00144 #endif
00145     return output ;
00146 }
00147 
00148 
00149 G4ITReactionChange* G4DNAMolecularReaction::MakeReaction(const G4Track& trackA, const G4Track& trackB)
00150 {
00151     fpChanges = new G4ITReactionChange();
00152     fpChanges->Initialize(trackA, trackB);
00153 
00154     G4Molecule* moleculeA = GetMolecule(trackA);
00155     G4Molecule* moleculeB = GetMolecule(trackB);
00156 
00157 #ifdef G4VERBOSE
00158     // DEBUG
00159     if(fVerbose)
00160     {
00161         G4cout << "G4DNAMolecularReaction::MakeReaction" << G4endl;
00162         G4cout<<"TrackA n°"                     << trackA.GetTrackID()
00163              <<"\t | Track B n°"                << trackB.GetTrackID() << G4endl;
00164 
00165         G4cout<<"Track A : Position : "         << G4BestUnit(trackA.GetPosition(),"Length")
00166              <<"\t Global Time : "              << G4BestUnit(trackA.GetGlobalTime(), "Time")<< G4endl;
00167 
00168         G4cout<<"Track B : Position : "         << G4BestUnit(trackB.GetPosition() ,"Length")
00169              <<"\t Global Time : "              << G4BestUnit(trackB.GetGlobalTime(), "Time")<< G4endl;
00170 
00171         G4cout<<"Reaction range : "             << G4BestUnit(fReactionRadius,"Length")
00172              << " \t Separation distance : "    << G4BestUnit(fDistance,"Length")<< G4endl;
00173         G4cout << "--------------------------------------------" << G4endl;
00174     }
00175 #endif
00176 
00177     const G4DNAMolecularReactionData* reactionData =  fMolReactionTable->GetReactionData(moleculeA, moleculeB);
00178 
00179     G4int nbProducts = reactionData->GetNbProducts();
00180 
00181     if (nbProducts)
00182     {
00183         G4double D1 = moleculeA->GetDiffusionCoefficient();
00184         G4double D2 = moleculeB->GetDiffusionCoefficient();
00185         G4double sqrD1 = sqrt(D1);
00186         G4double sqrD2 = sqrt(D2);
00187         G4double numerator = sqrD1 + sqrD2;
00188         G4ThreeVector reactionSite = sqrD1/numerator * trackA.GetPosition()
00189                                    + sqrD2/numerator * trackB.GetPosition()  ;
00190 
00191         for (G4int j=0 ; j < nbProducts; j++)
00192         {
00193             G4Molecule* product = new G4Molecule(*reactionData->GetProduct(j));
00194             G4Track* productTrack = product->BuildTrack(trackA.GetGlobalTime(),
00195                                                         reactionSite);
00196 
00197             productTrack->SetTrackStatus(fAlive);
00198 
00199             fpChanges->AddSecondary(productTrack);
00200             G4ITManager<G4Molecule>::Instance()->Push(productTrack);
00201         }
00202     }
00203 
00204     fpChanges->KillParents(true);
00205 
00206     return fpChanges;
00207 }

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