[1197] | 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 | // $Id: G4AdjointProcessEquivalentToDirectProcess.cc,v 1.1 2009/11/11 00:31:19 ldesorgh Exp $ |
---|
[1228] | 28 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[1197] | 29 | // |
---|
| 30 | // |
---|
| 31 | // ------------------------------------------------------------ |
---|
| 32 | // GEANT 4 class implementation file |
---|
| 33 | // |
---|
| 34 | // Class Description |
---|
| 35 | // |
---|
| 36 | // This class is for adjoint process equivalent to direct process |
---|
| 37 | |
---|
| 38 | // ------------------------------------------------------------ |
---|
| 39 | // Created by L.Desorgher 25 Sept. 2009 Inspired from G4WrapperProcess |
---|
| 40 | // ------------------------------------------------------------ |
---|
| 41 | |
---|
| 42 | #include "G4AdjointProcessEquivalentToDirectProcess.hh" |
---|
| 43 | #include "G4DynamicParticle.hh" |
---|
| 44 | G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess(const G4String& aName, |
---|
| 45 | G4VProcess* aProcess, |
---|
| 46 | G4ParticleDefinition* fwd_particle_def) |
---|
| 47 | :G4VProcess(aName) |
---|
| 48 | { |
---|
| 49 | theDirectProcess =aProcess; |
---|
| 50 | theProcessType = theDirectProcess->GetProcessType(); |
---|
| 51 | theFwdParticleDef = fwd_particle_def; |
---|
| 52 | } |
---|
| 53 | |
---|
| 54 | |
---|
| 55 | G4AdjointProcessEquivalentToDirectProcess::~G4AdjointProcessEquivalentToDirectProcess() |
---|
| 56 | { |
---|
| 57 | if (theDirectProcess!=0) delete theDirectProcess; |
---|
| 58 | } |
---|
| 59 | |
---|
| 60 | void G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft() |
---|
| 61 | { |
---|
| 62 | theDirectProcess->ResetNumberOfInteractionLengthLeft(); |
---|
| 63 | } |
---|
| 64 | |
---|
| 65 | G4double G4AdjointProcessEquivalentToDirectProcess:: |
---|
| 66 | AlongStepGetPhysicalInteractionLength( const G4Track& track, |
---|
| 67 | G4double previousStepSize, |
---|
| 68 | G4double currentMinimumStep, |
---|
| 69 | G4double& proposedSafety, |
---|
| 70 | G4GPILSelection* selection ) |
---|
| 71 | { |
---|
| 72 | |
---|
| 73 | |
---|
| 74 | //Change the particle definition to the direct one |
---|
| 75 | //------------------------------------------------ |
---|
| 76 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 77 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 78 | |
---|
| 79 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 80 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 81 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 82 | |
---|
| 83 | |
---|
| 84 | //Call the direct process |
---|
| 85 | //---------------------- |
---|
| 86 | G4double GPIL = theDirectProcess-> |
---|
| 87 | AlongStepGetPhysicalInteractionLength( track, |
---|
| 88 | previousStepSize, |
---|
| 89 | currentMinimumStep, |
---|
| 90 | proposedSafety, |
---|
| 91 | selection ); |
---|
| 92 | |
---|
| 93 | |
---|
| 94 | //Restore the adjoint particle definition to the direct one |
---|
| 95 | //------------------------------------------------ |
---|
| 96 | theDynPart->SetDefinition(adjPartDef); |
---|
| 97 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 98 | |
---|
| 99 | |
---|
| 100 | return GPIL; |
---|
| 101 | |
---|
| 102 | } |
---|
| 103 | |
---|
| 104 | G4double G4AdjointProcessEquivalentToDirectProcess:: |
---|
| 105 | AtRestGetPhysicalInteractionLength( const G4Track& track, |
---|
| 106 | G4ForceCondition* condition ) |
---|
| 107 | { //Change the particle definition to the direct one |
---|
| 108 | //------------------------------------------------ |
---|
| 109 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 110 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 111 | |
---|
| 112 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 113 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 114 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 115 | |
---|
| 116 | |
---|
| 117 | //Call the direct process |
---|
| 118 | //---------------------- |
---|
| 119 | |
---|
| 120 | |
---|
| 121 | G4double GPIL = theDirectProcess->AtRestGetPhysicalInteractionLength( track, condition ); |
---|
| 122 | |
---|
| 123 | //Restore the adjoint particle definition to the direct one |
---|
| 124 | //------------------------------------------------ |
---|
| 125 | theDynPart->SetDefinition(adjPartDef); |
---|
| 126 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 127 | |
---|
| 128 | return GPIL; |
---|
| 129 | |
---|
| 130 | |
---|
| 131 | } |
---|
| 132 | |
---|
| 133 | G4double G4AdjointProcessEquivalentToDirectProcess:: |
---|
| 134 | PostStepGetPhysicalInteractionLength( const G4Track& track, |
---|
| 135 | G4double previousStepSize, |
---|
| 136 | G4ForceCondition* condition ) |
---|
| 137 | { |
---|
| 138 | //Change the particle definition to the direct one |
---|
| 139 | //------------------------------------------------ |
---|
| 140 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 141 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 142 | |
---|
| 143 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 144 | |
---|
| 145 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 146 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 147 | |
---|
| 148 | |
---|
| 149 | //Call the direct process |
---|
| 150 | //---------------------- |
---|
| 151 | |
---|
| 152 | |
---|
| 153 | G4double GPIL = theDirectProcess->PostStepGetPhysicalInteractionLength( track, |
---|
| 154 | previousStepSize, |
---|
| 155 | condition ); |
---|
| 156 | |
---|
| 157 | //Restore the adjoint particle definition to the direct one |
---|
| 158 | //------------------------------------------------ |
---|
| 159 | theDynPart->SetDefinition(adjPartDef); |
---|
| 160 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 161 | |
---|
| 162 | return GPIL; |
---|
| 163 | |
---|
| 164 | |
---|
| 165 | } |
---|
| 166 | /* |
---|
| 167 | |
---|
| 168 | void G4AdjointProcessEquivalentToDirectProcess::SetProcessManager(const G4ProcessManager* procMan) |
---|
| 169 | { |
---|
| 170 | theDirectProcess->SetProcessManager(procMan); |
---|
| 171 | } |
---|
| 172 | |
---|
| 173 | const G4ProcessManager* G4AdjointProcessEquivalentToDirectProcess::GetProcessManager() |
---|
| 174 | { |
---|
| 175 | return theDirectProcess->GetProcessManager(); |
---|
| 176 | } |
---|
| 177 | */ |
---|
| 178 | G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt( const G4Track& track, |
---|
| 179 | const G4Step& stepData ) |
---|
| 180 | { |
---|
| 181 | //Change the particle definition to the direct one |
---|
| 182 | //------------------------------------------------ |
---|
| 183 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 184 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 185 | |
---|
| 186 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 187 | |
---|
| 188 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 189 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 190 | |
---|
| 191 | |
---|
| 192 | //Call the direct process |
---|
| 193 | //---------------------- |
---|
| 194 | |
---|
| 195 | G4VParticleChange* partChange = theDirectProcess->PostStepDoIt( track, stepData ); |
---|
| 196 | |
---|
| 197 | |
---|
| 198 | //Restore the adjoint particle definition to the direct one |
---|
| 199 | //------------------------------------------------ |
---|
| 200 | theDynPart->SetDefinition(adjPartDef); |
---|
| 201 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 202 | |
---|
| 203 | return partChange; |
---|
| 204 | |
---|
| 205 | |
---|
| 206 | |
---|
| 207 | } |
---|
| 208 | |
---|
| 209 | G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt( const G4Track& track, |
---|
| 210 | const G4Step& stepData ) |
---|
| 211 | { |
---|
| 212 | //Change the particle definition to the direct one |
---|
| 213 | //------------------------------------------------ |
---|
| 214 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 215 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 216 | |
---|
| 217 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 218 | |
---|
| 219 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 220 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 221 | |
---|
| 222 | |
---|
| 223 | //Call the direct process |
---|
| 224 | //---------------------- |
---|
| 225 | G4VParticleChange* partChange =theDirectProcess->AlongStepDoIt( track, stepData ); |
---|
| 226 | |
---|
| 227 | //Restore the adjoint particle definition to the direct one |
---|
| 228 | //------------------------------------------------ |
---|
| 229 | theDynPart->SetDefinition(adjPartDef); |
---|
| 230 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 231 | |
---|
| 232 | return partChange; |
---|
| 233 | } |
---|
| 234 | |
---|
| 235 | G4VParticleChange* G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt( const G4Track& track, |
---|
| 236 | const G4Step& stepData ) |
---|
| 237 | { |
---|
| 238 | //Change the particle definition to the direct one |
---|
| 239 | //------------------------------------------------ |
---|
| 240 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle()); |
---|
| 241 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 242 | |
---|
| 243 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 244 | |
---|
| 245 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 246 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 247 | |
---|
| 248 | |
---|
| 249 | //Call the direct process |
---|
| 250 | //---------------------- |
---|
| 251 | G4VParticleChange* partChange =theDirectProcess->AtRestDoIt( track, stepData ); |
---|
| 252 | |
---|
| 253 | //Restore the adjoint particle definition to the direct one |
---|
| 254 | //------------------------------------------------ |
---|
| 255 | theDynPart->SetDefinition(adjPartDef); |
---|
| 256 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 257 | |
---|
| 258 | return partChange; |
---|
| 259 | |
---|
| 260 | |
---|
| 261 | } |
---|
| 262 | |
---|
| 263 | G4bool G4AdjointProcessEquivalentToDirectProcess::IsApplicable(const G4ParticleDefinition&) |
---|
| 264 | { |
---|
| 265 | return theDirectProcess->IsApplicable(*theFwdParticleDef); |
---|
| 266 | } |
---|
| 267 | |
---|
| 268 | void G4AdjointProcessEquivalentToDirectProcess::BuildPhysicsTable(const G4ParticleDefinition& ) |
---|
| 269 | { |
---|
| 270 | return theDirectProcess->BuildPhysicsTable(*theFwdParticleDef); |
---|
| 271 | } |
---|
| 272 | |
---|
| 273 | void G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable(const G4ParticleDefinition& ) |
---|
| 274 | { |
---|
| 275 | return theDirectProcess->PreparePhysicsTable(*theFwdParticleDef); |
---|
| 276 | } |
---|
| 277 | |
---|
| 278 | G4bool G4AdjointProcessEquivalentToDirectProcess:: |
---|
| 279 | StorePhysicsTable(const G4ParticleDefinition* , |
---|
| 280 | const G4String& directory, |
---|
| 281 | G4bool ascii) |
---|
| 282 | { |
---|
| 283 | return theDirectProcess->StorePhysicsTable(theFwdParticleDef, directory, ascii); |
---|
| 284 | } |
---|
| 285 | |
---|
| 286 | G4bool G4AdjointProcessEquivalentToDirectProcess:: |
---|
| 287 | RetrievePhysicsTable( const G4ParticleDefinition* , |
---|
| 288 | const G4String& directory, |
---|
| 289 | G4bool ascii) |
---|
| 290 | { |
---|
| 291 | return theDirectProcess->RetrievePhysicsTable(theFwdParticleDef, directory, ascii); |
---|
| 292 | } |
---|
| 293 | |
---|
| 294 | void G4AdjointProcessEquivalentToDirectProcess::StartTracking(G4Track* track) |
---|
| 295 | { |
---|
| 296 | //Change the particle definition to the direct one |
---|
| 297 | //------------------------------------------------ |
---|
| 298 | G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track->GetDynamicParticle()); |
---|
| 299 | G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition(); |
---|
| 300 | |
---|
| 301 | G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts()); |
---|
| 302 | theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0)); |
---|
| 303 | theDynPart->SetDefinition(theFwdParticleDef); |
---|
| 304 | |
---|
| 305 | theDirectProcess->StartTracking(track); |
---|
| 306 | |
---|
| 307 | //Restore the adjoint particle definition to the direct one |
---|
| 308 | //------------------------------------------------ |
---|
| 309 | theDynPart->SetDefinition(adjPartDef); |
---|
| 310 | theDynPart->SetPreAssignedDecayProducts(decayProducts); |
---|
| 311 | |
---|
| 312 | |
---|
| 313 | return; |
---|
| 314 | |
---|
| 315 | } |
---|
| 316 | |
---|
| 317 | void G4AdjointProcessEquivalentToDirectProcess::EndTracking() |
---|
| 318 | { |
---|
| 319 | theDirectProcess->EndTracking(); |
---|
| 320 | } |
---|
| 321 | |
---|