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 #include "G4DNAUpdateSystemModel.hh" 28 #include "G4Molecule.hh" 29 #include "G4DNAMolecularReactionTable.hh" 30 #include "G4UnitsTable.hh" 31 #include "G4MoleculeCounter.hh" 32 #include "G4DNAScavengerMaterial.hh" 33 #include "G4Scheduler.hh" 34 35 G4DNAUpdateSystemModel::G4DNAUpdateSystemModel() = default; 36 37 void G4DNAUpdateSystemModel::SetMesh(G4DNAMesh* pMesh) { fpMesh = pMesh; } 38 void G4DNAUpdateSystemModel::KillMolecule(const Index& index, MolType type) 39 { 40 // kill normal molecule 41 auto& node = fpMesh->GetVoxelMapList(index); 42 auto iter = node.find(type); 43 if(iter != node.end()) 44 { 45 if(iter->second <= 0) 46 { 47 G4ExceptionDescription exceptionDescription; 48 exceptionDescription 49 << "G4DNAUpdateSystemModel::KillMolecule::molecule : " 50 << type->GetName() << " index : " << index 51 << " number : " << iter->second << G4endl; 52 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler002", 53 FatalErrorInArgument, exceptionDescription); 54 } 55 iter->second--; 56 if(G4VMoleculeCounter::Instance()->InUse()) 57 { 58 G4VMoleculeCounter::Instance()->RemoveAMoleculeAtTime(type, fGlobalTime); 59 } 60 } 61 else 62 { 63 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>( 64 G4Scheduler::Instance()->GetScavengerMaterial()); 65 if(pScavengerMaterial != nullptr) 66 { 67 pScavengerMaterial->ReduceNumberMoleculePerVolumeUnitForMaterialConf( 68 type, fGlobalTime); 69 } 70 else 71 { 72 G4ExceptionDescription exceptionDescription; 73 exceptionDescription 74 << "index : " << index << " " << type->GetName() 75 << " This molecule is not belong scavengers or particle-base" 76 << G4endl; 77 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler002", 78 FatalErrorInArgument, exceptionDescription); 79 } 80 } 81 } 82 83 void G4DNAUpdateSystemModel::JumpTo(const Index& index, MolType type) 84 { 85 auto& node = fpMesh->GetVoxelMapList(index); 86 auto iter = node.find(type); 87 if(iter != node.end()) 88 { 89 if(iter->second <= 0) 90 { 91 G4ExceptionDescription exceptionDescription; 92 exceptionDescription << "G4DNAUpdateSystemModel::JumpTo::molecule : " 93 << type->GetName() << " index : " << index 94 << " number : " << iter->second; 95 G4Exception("G4DNAUpdateSystemModel::JumpTo", "G4DNAUpdateSystemModel001", 96 FatalErrorInArgument, exceptionDescription); 97 } 98 iter->second--; 99 } 100 else 101 { 102 fpMesh->PrintVoxel(index); 103 G4ExceptionDescription exceptionDescription; 104 exceptionDescription << "index : " << index << " " << type->GetName() 105 << " There is no this type"; 106 G4Exception("G4DNAUpdateSystemModel::JumpTo", "G4DNAUpdateSystemModel002", 107 FatalErrorInArgument, exceptionDescription); 108 } 109 } 110 111 void G4DNAUpdateSystemModel::CreateMolecule(const Index& index, MolType type) 112 { 113 // for scavenger 114 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>( 115 G4Scheduler::Instance()->GetScavengerMaterial()); 116 if(pScavengerMaterial != nullptr && pScavengerMaterial->find(type)) 117 { 118 pScavengerMaterial->AddNumberMoleculePerVolumeUnitForMaterialConf( 119 type, fGlobalTime); 120 return; 121 } 122 // for molecule 123 auto& node = fpMesh->GetVoxelMapList(index); 124 auto iter = node.find(type); 125 if(iter != node.end()) 126 { 127 iter->second++; 128 } 129 else 130 { 131 node[type] = 1; 132 } 133 134 if(G4VMoleculeCounter::Instance()->InUse()) 135 { 136 G4VMoleculeCounter::Instance()->AddAMoleculeAtTime(type, fGlobalTime); 137 } 138 } 139 140 void G4DNAUpdateSystemModel::JumpIn(const Index& index, MolType type) 141 { 142 // for molecule 143 auto& node = fpMesh->GetVoxelMapList(index); 144 auto iter = node.find(type); 145 if(iter != node.end()) 146 { 147 iter->second++; 148 } 149 else 150 { 151 node[type] = 1; 152 } 153 } 154 155 void G4DNAUpdateSystemModel::UpdateSystem(const Index& index, 156 const ReactionData& data) 157 { 158 auto reactant1 = data.GetReactant1(); 159 auto reactant2 = data.GetReactant2(); 160 #ifdef G4VERBOSE 161 if(fVerbose != 0) 162 { 163 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time") 164 << " Reaction : " << reactant1->GetName() << " + " 165 << reactant2->GetName() << " -> "; 166 } 167 #endif 168 const G4int nbProducts = data.GetNbProducts(); 169 if(nbProducts != 0) 170 { 171 for(G4int j = 0; j < nbProducts; ++j) 172 { 173 #ifdef G4VERBOSE 174 if((fVerbose != 0) && j != 0) 175 { 176 G4cout << " + "; 177 } 178 if(fVerbose != 0) 179 { 180 G4cout << data.GetProduct(j)->GetName(); 181 } 182 #endif 183 CreateMolecule(index, data.GetProduct(j)); 184 } 185 } 186 else 187 { 188 #ifdef G4VERBOSE 189 if(fVerbose != 0) 190 { 191 G4cout << "No product"; 192 } 193 #endif 194 } 195 #ifdef G4VERBOSE 196 if(fVerbose != 0) 197 { 198 G4cout << G4endl; 199 } 200 #endif 201 KillMolecule(index, reactant1); 202 KillMolecule(index, reactant2); 203 } 204 205 void G4DNAUpdateSystemModel::UpdateSystem(const Index& index, 206 const JumpingData& data) 207 { 208 auto reactant = std::get<0>(data); 209 auto JunpToIndex = std::get<1>(data); 210 #ifdef G4VERBOSE 211 if(fVerbose > 1) 212 { 213 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time") 214 << " Jumping : " << reactant->GetName() << " from " << index 215 << " -> " << JunpToIndex << G4endl; 216 } 217 #endif 218 JumpTo(index, reactant); 219 JumpIn(JunpToIndex, reactant); 220 } 221