G4AdjointProcessEquivalentToDirectProcess.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 //
00027 // $Id: G4AdjointProcessEquivalentToDirectProcess.cc 69844 2013-05-16 09:19:33Z gcosmo $
00028 //
00029 // 
00030 // ------------------------------------------------------------
00031 //        GEANT 4 class implementation file 
00032 //
00033 // Class Description
00034 //
00035 // This class is for adjoint process equivalent to direct process
00036 
00037 // ------------------------------------------------------------
00038 //   Created by L.Desorgher          25 Sept. 2009  Inspired from G4WrapperProcess
00039 // ------------------------------------------------------------
00040 
00041 #include "G4AdjointProcessEquivalentToDirectProcess.hh"
00042 #include "G4DynamicParticle.hh"
00043 G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess(const G4String& aName,
00044                                                                                      G4VProcess* aProcess,
00045                                                                                      G4ParticleDefinition* fwd_particle_def)
00046 :G4VProcess(aName)
00047 {  
00048    theDirectProcess =aProcess;
00049    theProcessType = theDirectProcess->GetProcessType();
00050    theFwdParticleDef  = fwd_particle_def;
00051 }
00052 
00053 
00054 G4AdjointProcessEquivalentToDirectProcess::~G4AdjointProcessEquivalentToDirectProcess()
00055 {
00056   if (theDirectProcess!=0) delete theDirectProcess;
00057 }
00058 
00059 void G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft()
00060 {
00061   theDirectProcess->ResetNumberOfInteractionLengthLeft();
00062 }
00063 
00064 G4double G4AdjointProcessEquivalentToDirectProcess::
00065 AlongStepGetPhysicalInteractionLength( const G4Track& track,
00066                                              G4double  previousStepSize,
00067                                              G4double  currentMinimumStep,
00068                                              G4double& proposedSafety,
00069                                              G4GPILSelection* selection     )
00070 { 
00071   
00072 
00073   //Change the particle definition to the direct one
00074   //------------------------------------------------
00075   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00076   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00077   
00078   G4DecayProducts* decayProducts = const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00079   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00080   theDynPart->SetDefinition(theFwdParticleDef);
00081   
00082   
00083   //Call the direct process
00084   //----------------------
00085   G4double GPIL =  theDirectProcess->
00086          AlongStepGetPhysicalInteractionLength( track,
00087                                                 previousStepSize,
00088                                                 currentMinimumStep,
00089                                                 proposedSafety,
00090                                                 selection     );
00091   
00092   
00093   //Restore the adjoint particle definition to the direct one
00094   //------------------------------------------------
00095   theDynPart->SetDefinition(adjPartDef);
00096   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00097   
00098   
00099   return GPIL;
00100                                                 
00101 }
00102 
00103 G4double G4AdjointProcessEquivalentToDirectProcess::
00104 AtRestGetPhysicalInteractionLength( const G4Track& track,
00105                                           G4ForceCondition* condition )
00106 { //Change the particle definition to the direct one
00107   //------------------------------------------------
00108   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00109   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00110   
00111   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00112   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00113   theDynPart->SetDefinition(theFwdParticleDef);
00114   
00115   
00116   //Call the direct process
00117   //----------------------
00118 
00119    
00120   G4double GPIL =  theDirectProcess->AtRestGetPhysicalInteractionLength( track, condition );
00121   
00122   //Restore the adjoint particle definition to the direct one
00123   //------------------------------------------------
00124   theDynPart->SetDefinition(adjPartDef);
00125   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00126   
00127   return GPIL;
00128                                                 
00129   
00130 }
00131 
00132 G4double G4AdjointProcessEquivalentToDirectProcess::
00133 PostStepGetPhysicalInteractionLength( const G4Track& track,
00134                                             G4double   previousStepSize,
00135                                             G4ForceCondition* condition )
00136 {
00137   //Change the particle definition to the direct one
00138   //------------------------------------------------
00139   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00140   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00141   
00142   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00143  
00144   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00145   theDynPart->SetDefinition(theFwdParticleDef);
00146   
00147   
00148   //Call the direct process
00149   //----------------------
00150 
00151    
00152   G4double GPIL = theDirectProcess->PostStepGetPhysicalInteractionLength( track,
00153                                                              previousStepSize,
00154                                                              condition );
00155    
00156   //Restore the adjoint particle definition to the direct one
00157   //------------------------------------------------
00158   theDynPart->SetDefinition(adjPartDef);
00159   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00160   
00161    return GPIL;
00162   
00163                                      
00164 }
00165 /*
00166       
00167 void G4AdjointProcessEquivalentToDirectProcess::SetProcessManager(const G4ProcessManager* procMan)
00168 {
00169    theDirectProcess->SetProcessManager(procMan); 
00170 }
00171 
00172 const G4ProcessManager* G4AdjointProcessEquivalentToDirectProcess::GetProcessManager()
00173 {
00174   return     theDirectProcess->GetProcessManager();
00175 }
00176 */
00177 G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt( const G4Track& track,
00178                                                    const G4Step&  stepData )
00179 {
00180   //Change the particle definition to the direct one
00181   //------------------------------------------------
00182   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00183   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00184   
00185   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00186  
00187   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00188   theDynPart->SetDefinition(theFwdParticleDef);
00189   
00190   
00191   //Call the direct process
00192   //----------------------
00193   
00194   G4VParticleChange* partChange = theDirectProcess->PostStepDoIt( track, stepData );
00195   
00196   
00197   //Restore the adjoint particle definition to the direct one
00198   //------------------------------------------------
00199   theDynPart->SetDefinition(adjPartDef);
00200   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00201   
00202   return partChange;
00203   
00204   
00205           
00206 }
00207 
00208 G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt( const G4Track& track,
00209                                                     const G4Step& stepData )
00210 { 
00211   //Change the particle definition to the direct one
00212   //------------------------------------------------
00213   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00214   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00215   
00216   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00217  
00218   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00219   theDynPart->SetDefinition(theFwdParticleDef);
00220   
00221   
00222   //Call the direct process
00223   //----------------------
00224   G4VParticleChange* partChange =theDirectProcess->AlongStepDoIt( track, stepData );
00225   
00226   //Restore the adjoint particle definition to the direct one
00227   //------------------------------------------------
00228   theDynPart->SetDefinition(adjPartDef);
00229   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00230   
00231   return partChange;        
00232 }
00233  
00234 G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt( const G4Track& track,
00235                                                  const G4Step& stepData )
00236 {
00237   //Change the particle definition to the direct one
00238   //------------------------------------------------
00239   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
00240   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00241   
00242   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*>  (theDynPart->GetPreAssignedDecayProducts());
00243  
00244   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00245   theDynPart->SetDefinition(theFwdParticleDef);
00246   
00247   
00248   //Call the direct process
00249   //----------------------
00250   G4VParticleChange* partChange =theDirectProcess->AtRestDoIt( track, stepData );
00251   
00252   //Restore the adjoint particle definition to the direct one
00253   //------------------------------------------------
00254   theDynPart->SetDefinition(adjPartDef);
00255   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00256   
00257    return partChange; 
00258   
00259        
00260 }
00261 
00262 G4bool G4AdjointProcessEquivalentToDirectProcess::IsApplicable(const G4ParticleDefinition&)
00263 {
00264   return     theDirectProcess->IsApplicable(*theFwdParticleDef);
00265 }
00266 
00267 void G4AdjointProcessEquivalentToDirectProcess::BuildPhysicsTable(const G4ParticleDefinition& )
00268 {
00269   return     theDirectProcess->BuildPhysicsTable(*theFwdParticleDef);
00270 }
00271 
00272 void G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable(const G4ParticleDefinition& )
00273 {
00274   return     theDirectProcess->PreparePhysicsTable(*theFwdParticleDef);
00275 }
00276 
00277 G4bool G4AdjointProcessEquivalentToDirectProcess::
00278 StorePhysicsTable(const G4ParticleDefinition* ,
00279                   const G4String& directory, 
00280                         G4bool          ascii)
00281 {
00282   return theDirectProcess->StorePhysicsTable(theFwdParticleDef,  directory,  ascii);
00283 } 
00284  
00285 G4bool G4AdjointProcessEquivalentToDirectProcess::
00286 RetrievePhysicsTable( const G4ParticleDefinition* ,
00287                       const G4String& directory, 
00288                             G4bool          ascii)
00289 {
00290   return theDirectProcess->RetrievePhysicsTable(theFwdParticleDef,  directory,  ascii);
00291 }  
00292 
00293 void G4AdjointProcessEquivalentToDirectProcess::StartTracking(G4Track* track)
00294 {
00295   //Change the particle definition to the direct one
00296   //------------------------------------------------
00297   G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track->GetDynamicParticle());
00298   G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
00299   
00300   G4DecayProducts* decayProducts =  const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
00301   theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
00302   theDynPart->SetDefinition(theFwdParticleDef);
00303   
00304   theDirectProcess->StartTracking(track);
00305   
00306    //Restore the adjoint particle definition to the direct one
00307   //------------------------------------------------
00308   theDynPart->SetDefinition(adjPartDef);
00309   theDynPart->SetPreAssignedDecayProducts(decayProducts);
00310   
00311   
00312   return; 
00313   
00314 }
00315 
00316 void G4AdjointProcessEquivalentToDirectProcess::EndTracking()
00317 {
00318   theDirectProcess->EndTracking(); 
00319 }
00320 

Generated on Mon May 27 17:47:37 2013 for Geant4 by  doxygen 1.4.7