G4DNAMolecularReaction Class Reference

#include <G4DNAMolecularReaction.hh>

Inheritance diagram for G4DNAMolecularReaction:

G4VITReactionProcess

Public Member Functions

 G4DNAMolecularReaction ()
virtual ~G4DNAMolecularReaction ()
 G4DNAMolecularReaction (const G4DNAMolecularReaction &other)
G4DNAMolecularReactionoperator= (const G4DNAMolecularReaction &other)
virtual G4bool TestReactibility (const G4Track &, const G4Track &, const double currentStepTime, const double previousStepTime, bool userStepTimeLimit)
virtual G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &)
void SetReactionModel (G4VDNAReactionModel *)
void SetReactionTable (const G4DNAMolecularReactionTable *)
void SetVerbose (int)

Protected Attributes

const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModelfReactionModel
G4int fVerbose
G4double fReactionRadius
G4double fDistance

Detailed Description

G4DNAMolecularReaction is the reaction process used in G4DNAMolecularStepByStepModel between two molecules. After the global track steps, it test whether the molecules can react. If so, the reaction is made.

Definition at line 55 of file G4DNAMolecularReaction.hh.


Constructor & Destructor Documentation

G4DNAMolecularReaction::G4DNAMolecularReaction (  ) 

Default constructor

Definition at line 47 of file G4DNAMolecularReaction.cc.

References fDistance, fReactionModel, fReactionRadius, and fVerbose.

00047                                               :G4VITReactionProcess(),
00048     fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
00049 {
00050     //ctor
00051     fVerbose        = 0;
00052     fReactionModel  = 0;
00053     fReactionRadius = -1;
00054     fDistance       = -1;
00055 }

G4DNAMolecularReaction::~G4DNAMolecularReaction (  )  [virtual]

Default destructor

Definition at line 57 of file G4DNAMolecularReaction.cc.

References G4VITReactionProcess::fpChanges.

00058 {
00059     //dtor
00060     if(fpChanges) delete fpChanges;
00061 }

G4DNAMolecularReaction::G4DNAMolecularReaction ( const G4DNAMolecularReaction other  ) 

Copy constructor

Parameters:
other Object to copy from

Definition at line 63 of file G4DNAMolecularReaction.cc.

References fDistance, fMolReactionTable, fReactionModel, fReactionRadius, and fVerbose.

00063                                                                                  :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 }


Member Function Documentation

G4ITReactionChange * G4DNAMolecularReaction::MakeReaction ( const G4Track ,
const G4Track  
) [virtual]

Will generate the products of the two given tracks

Implements G4VITReactionProcess.

Definition at line 149 of file G4DNAMolecularReaction.cc.

References G4ITReactionChange::AddSecondary(), fAlive, fDistance, fMolReactionTable, G4VITReactionProcess::fpChanges, fReactionRadius, fVerbose, G4BestUnit, G4cout, G4endl, G4Molecule::GetDiffusionCoefficient(), G4Track::GetGlobalTime(), GetMolecule(), G4DNAMolecularReactionData::GetNbProducts(), G4Track::GetPosition(), G4DNAMolecularReactionData::GetProduct(), G4DNAMolecularReactionTable::GetReactionData(), G4Track::GetTrackID(), G4ITReactionChange::Initialize(), G4ITManager< T >::Instance(), G4ITReactionChange::KillParents(), and G4Track::SetTrackStatus().

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 }

G4DNAMolecularReaction & G4DNAMolecularReaction::operator= ( const G4DNAMolecularReaction other  ) 

Assignment operator

Parameters:
other Object to assign from
Returns:
A reference to this

Definition at line 74 of file G4DNAMolecularReaction.cc.

References fDistance, fMolReactionTable, fReactionRadius, and fVerbose.

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 }

void G4DNAMolecularReaction::SetReactionModel ( G4VDNAReactionModel  )  [inline]

Definition at line 102 of file G4DNAMolecularReaction.hh.

References fReactionModel.

00103 {
00104     fReactionModel = model;
00105 }

void G4DNAMolecularReaction::SetReactionTable ( const G4DNAMolecularReactionTable  )  [inline]

void G4DNAMolecularReaction::SetVerbose ( int   )  [inline]

Definition at line 107 of file G4DNAMolecularReaction.hh.

References fVerbose.

00108 {
00109     fVerbose = verb;
00110 }

G4bool G4DNAMolecularReaction::TestReactibility ( const G4Track ,
const G4Track ,
const double  currentStepTime,
const double  previousStepTime,
bool  userStepTimeLimit 
) [virtual]

Given two tracks, it tells you whether they can react

Implements G4VITReactionProcess.

Definition at line 87 of file G4DNAMolecularReaction.cc.

References FatalErrorInArgument, fDistance, G4VDNAReactionModel::FindReaction(), fMolReactionTable, fReactionModel, fReactionRadius, fVerbose, G4BestUnit, G4cout, G4endl, G4Exception(), GetMolecule(), G4Molecule::GetName(), G4Track::GetPosition(), and G4Track::GetTrackID().

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 }


Field Documentation

G4double G4DNAMolecularReaction::fDistance [protected]

Definition at line 98 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), MakeReaction(), operator=(), and TestReactibility().

const G4DNAMolecularReactionTable*& G4DNAMolecularReaction::fMolReactionTable [protected]

Definition at line 94 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), MakeReaction(), operator=(), and TestReactibility().

G4VDNAReactionModel* G4DNAMolecularReaction::fReactionModel [protected]

Definition at line 95 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), SetReactionModel(), and TestReactibility().

G4double G4DNAMolecularReaction::fReactionRadius [protected]

Definition at line 97 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), MakeReaction(), operator=(), and TestReactibility().

G4int G4DNAMolecularReaction::fVerbose [protected]

Definition at line 96 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), MakeReaction(), operator=(), SetVerbose(), and TestReactibility().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:48 2013 for Geant4 by  doxygen 1.4.7