Geant4 Cross Reference |
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 /// \file PhysChemIO.cc 28 /// \brief Implementation of the PhysChemIO class 29 30 #include "PhysChemIO.hh" 31 #include "SteppingAction.hh" 32 #include "Analysis.hh" 33 #include "G4Track.hh" 34 #include "G4NavigationHistory.hh" 35 #include "G4RunManager.hh" 36 37 #ifdef USE_MPI 38 #include "G4MPImanager.hh" 39 #endif 40 41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 42 43 PhysChemIO::PhysChemIO(SteppingAction* stepAction) : G4VPhysChemIO(), 44 fSteppingAction(stepAction) 45 {;} 46 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 48 49 void PhysChemIO::CreateWaterMolecule(G4int electronicModif, G4int electronicLevel, 50 G4double /*energy*/, 51 const G4Track* theIncomingTrack) 52 { 53 //L.T. Anh: to correct electronicLevel in G4DNAChemistryManager, 54 //see in G4DNAChemistryManager::CreateWaterMolecule 55 electronicLevel = 4 - electronicLevel; 56 //Rel pos 57 G4ThreeVector relPos; 58 auto touchable = theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable(); 59 relPos = touchable->GetHistory()->GetTopTransform().TransformPoint(theIncomingTrack->GetPosition()); 60 61 // Get the flag of the current volume 62 G4int volumeFlag =fSteppingAction->SetupVolumeFlag( 63 theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()); 64 65 if( volumeFlag == 161 // voxelStraight 66 || volumeFlag == 162 // voxelRight 67 || volumeFlag == 163 // voxelLeft 68 || volumeFlag == 164 // voxelUp 69 || volumeFlag == 165 // voxelDown 70 || volumeFlag == 261 // voxelStraight2 71 || volumeFlag == 262 // voxelRight2 72 || volumeFlag == 263 // voxelLeft2 73 || volumeFlag == 264 // voxelUp2 74 || volumeFlag == 265) // voxelDown2 75 { 76 // Get the volume copy number 77 G4int volumeCpNum = touchable->GetCopyNumber(); 78 //theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetVolume()->GetUserID(); 79 80 // Default flag values 81 G4String motherVolumeName = ""; 82 G4int motherVolumeFlag = -1; 83 G4int motherVolumeCpNum = -1; 84 85 // Mother volume informations 86 87 // Be sure there is a mother volume to ask for 88 if(theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetHistoryDepth() >0) 89 { 90 G4VPhysicalVolume* motherVol = theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetVolume(1); 91 92 // General infos 93 motherVolumeName = motherVol->GetName(); 94 motherVolumeFlag = fSteppingAction->SetupVolumeFlag(motherVolumeName); 95 motherVolumeCpNum = motherVol->GetCopyNo(); 96 } 97 G4int eventId = G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID(); 98 #ifdef USE_MPI 99 auto g4MPI = G4MPImanager::GetManager(); 100 if (g4MPI->IsSlave()) { // update eventID only for slave, cause rank_master=0 101 G4int rank = g4MPI->GetRank(); 102 eventId += g4MPI->GetEventsInMaster() + (rank-1)*g4MPI->GetEventsInSlave(); 103 } 104 #endif 105 106 InfoForChemGeo aInfo; 107 aInfo.fType = 1; // water=1 108 aInfo.fState = electronicModif ; 109 aInfo.fElectronicLevel = electronicLevel ; 110 aInfo.fX = theIncomingTrack->GetPosition().x()/nm; 111 aInfo.fY = theIncomingTrack->GetPosition().y()/nm; 112 aInfo.fZ = theIncomingTrack->GetPosition().z()/nm; 113 aInfo.fParentTrackID = theIncomingTrack->GetTrackID() ; 114 aInfo.fEventNumber = eventId; 115 aInfo.fVolume = volumeFlag ; 116 aInfo.fVolumeCopyNumber = volumeCpNum; 117 aInfo.fMotherVolume = motherVolumeFlag ; 118 aInfo.fMotherVolumeCopyNumber = motherVolumeCpNum ; 119 aInfo.fRelX = relPos.x()/nm; 120 aInfo.fRelY = relPos.y()/nm; 121 aInfo.fRelZ = relPos.z()/nm; 122 Analysis::GetAnalysis()->AddInfoForChemGeo(aInfo); 123 } 124 } 125 126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 127 128 void PhysChemIO::CreateSolvatedElectron(const G4Track* theIncomingTrack, G4ThreeVector* finalPosition) 129 { 130 G4ThreeVector pos; 131 if(finalPosition) pos = *finalPosition; 132 else pos = theIncomingTrack->GetPosition(); 133 134 // Rel pos 135 G4ThreeVector relPos; 136 const G4VTouchable* touchable = theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable(); 137 relPos = touchable->GetHistory()->GetTopTransform().TransformPoint(pos); 138 139 // Current volume infos 140 G4int volumeFlag = fSteppingAction->SetupVolumeFlag( 141 theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()); 142 if( volumeFlag == 161 // voxelStraight 143 || volumeFlag == 162 // voxelRight 144 || volumeFlag == 163 // voxelLeft 145 || volumeFlag == 164 // voxelUp 146 || volumeFlag == 165 // voxelDown 147 || volumeFlag == 261 // voxelStraight2 148 || volumeFlag == 262 // voxelRight2 149 || volumeFlag == 263 // voxelLeft2 150 || volumeFlag == 264 // voxelUp2 151 || volumeFlag == 265) // voxelDown2 152 { 153 G4int volumeCpNum = touchable->GetCopyNumber(); 154 155 G4String motherVolumeName = ""; 156 G4int motherVolumeFlag = -1; 157 G4int motherVolumeCpNum = -1; 158 159 // Mother volume informations 160 // Be sure there is a mother volume to ask for 161 if(theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetHistoryDepth() >0) 162 { 163 G4VPhysicalVolume* motherVol = theIncomingTrack->GetStep()->GetPreStepPoint()->GetTouchable()->GetVolume(1); 164 // General infos 165 motherVolumeName = motherVol->GetName(); 166 motherVolumeFlag = fSteppingAction->SetupVolumeFlag(motherVolumeName); 167 motherVolumeCpNum = motherVol->GetCopyNo(); 168 } 169 G4int eventId = G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID(); 170 #ifdef USE_MPI 171 auto g4MPI = G4MPImanager::GetManager(); 172 if (g4MPI->IsSlave()) { // update eventID only for slave, cause rank_master=0 173 G4int rank = g4MPI->GetRank(); 174 eventId += g4MPI->GetEventsInMaster() + (rank-1)*g4MPI->GetEventsInSlave(); 175 } 176 #endif 177 178 InfoForChemGeo aInfo; 179 aInfo.fType = 2; // / solvated electron=2 180 aInfo.fState = -1; // no state for solvated electron 181 aInfo.fElectronicLevel = -1; // no electronic level for solvated electron 182 aInfo.fX = pos.x()/nm; 183 aInfo.fY = pos.y()/nm; 184 aInfo.fZ = pos.z()/nm; 185 aInfo.fParentTrackID = theIncomingTrack->GetTrackID() ; 186 aInfo.fEventNumber = eventId; 187 aInfo.fVolume = volumeFlag ; 188 aInfo.fVolumeCopyNumber = volumeCpNum ; 189 aInfo.fMotherVolume = motherVolumeFlag ; 190 aInfo.fMotherVolumeCopyNumber = motherVolumeCpNum ; 191 aInfo.fRelX = relPos.x()/nm; 192 aInfo.fRelY = relPos.y()/nm; 193 aInfo.fRelZ = relPos.z()/nm; 194 Analysis::GetAnalysis()->AddInfoForChemGeo(aInfo); 195 } 196 } 197 198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......