[819] | 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 | #include "G4QMDCollision.hh" |
---|
| 27 | #include "G4ParticleDefinition.hh" |
---|
| 28 | #include "G4Scatterer.hh" |
---|
| 29 | #include "Randomize.hh" |
---|
| 30 | #include "G4PionZero.hh" |
---|
| 31 | |
---|
| 32 | G4QMDCollision::G4QMDCollision() |
---|
| 33 | : deltar ( 4 ) |
---|
| 34 | , bcmax0 ( 1.323142 ) // NN maximum impact parameter |
---|
| 35 | , bcmax1 ( 2.523 ) // others maximum impact parameter |
---|
| 36 | , sig0 ( 55 ) // NN cross section |
---|
| 37 | , sig1 ( 200 ) // others cross section |
---|
| 38 | , epse ( 0.0001 ) |
---|
| 39 | { |
---|
| 40 | theScatterer = new G4Scatterer(); |
---|
| 41 | } |
---|
| 42 | |
---|
| 43 | |
---|
| 44 | |
---|
| 45 | G4QMDCollision::~G4QMDCollision() |
---|
| 46 | { |
---|
| 47 | delete theScatterer; |
---|
| 48 | } |
---|
| 49 | |
---|
| 50 | |
---|
| 51 | void G4QMDCollision::CalKinematicsOfBinaryCollisions() |
---|
| 52 | { |
---|
| 53 | |
---|
| 54 | |
---|
| 55 | G4int n = theSystem->GetTotalNumberOfParticipant(); |
---|
| 56 | |
---|
| 57 | //071101 |
---|
| 58 | for ( G4int i = 0 ; i < n ; i++ ) |
---|
| 59 | { |
---|
| 60 | //std::cout << i << " " << theSystem->GetParticipant( i )->GetDefinition()->GetParticleName() << " " << theSystem->GetParticipant( i )->GetPosition() << std::endl; |
---|
| 61 | if ( theSystem->GetParticipant( i )->GetDefinition()->IsShortLived() ) |
---|
| 62 | { |
---|
| 63 | G4ParticleDefinition* pd0 = theSystem->GetParticipant( i )->GetDefinition(); |
---|
| 64 | G4ThreeVector p0 = theSystem->GetParticipant( i )->GetMomentum(); |
---|
| 65 | G4ThreeVector r0 = theSystem->GetParticipant( i )->GetPosition(); |
---|
| 66 | |
---|
| 67 | G4LorentzVector p40 = theSystem->GetParticipant( i )->Get4Momentum(); |
---|
| 68 | |
---|
| 69 | G4double epot = theMeanField->GetTotalPotential(); |
---|
| 70 | G4double eini = epot + p40.e(); |
---|
| 71 | |
---|
| 72 | G4int n0 = theSystem->GetTotalNumberOfParticipant(); |
---|
| 73 | G4int i0 = 0; |
---|
| 74 | G4bool isThisEnergyOK = false; |
---|
| 75 | for ( G4int ii = 0 ; ii < 4 ; ii++ ) |
---|
| 76 | { |
---|
| 77 | |
---|
| 78 | //G4LorentzVector p4 = theSystem->GetParticipant( i )->Get4Momentum(); |
---|
| 79 | G4LorentzVector p400 = p40; |
---|
| 80 | |
---|
| 81 | p400 *= GeV; |
---|
| 82 | //G4KineticTrack kt( theSystem->GetParticipant( i )->GetDefinition() , 0.0 , (theSystem->GetParticipant( i )->GetPosition())*fermi , p4 ); |
---|
| 83 | G4KineticTrack kt( pd0 , 0.0 , r0*fermi , p400 ); |
---|
| 84 | // std::cout << "G4KineticTrack " << i << " " << kt.GetDefinition()->GetParticleName() << kt.GetPosition() << std::endl; |
---|
| 85 | G4KineticTrackVector* secs = NULL; |
---|
| 86 | secs = kt.Decay(); |
---|
| 87 | G4int id = 0; |
---|
| 88 | G4double et = 0; |
---|
| 89 | if ( secs ) |
---|
| 90 | { |
---|
| 91 | for ( G4KineticTrackVector::iterator it |
---|
| 92 | = secs->begin() ; it != secs->end() ; it++ ) |
---|
| 93 | { |
---|
| 94 | // std::cout << "G4KineticTrack" |
---|
| 95 | // << " " << (*it)->GetDefinition()->GetParticleName() |
---|
| 96 | // << " " << (*it)->Get4Momentum() |
---|
| 97 | // << " " << (*it)->GetPosition()/fermi |
---|
| 98 | // << std::endl; |
---|
| 99 | if ( id == 0 ) |
---|
| 100 | { |
---|
| 101 | theSystem->GetParticipant( i )->SetDefinition( (*it)->GetDefinition() ); |
---|
| 102 | theSystem->GetParticipant( i )->SetMomentum( (*it)->Get4Momentum().v()/GeV ); |
---|
| 103 | theSystem->GetParticipant( i )->SetPosition( (*it)->GetPosition()/fermi ); |
---|
| 104 | //theMeanField->Cal2BodyQuantities( i ); |
---|
| 105 | et += (*it)->Get4Momentum().e()/GeV; |
---|
| 106 | } |
---|
| 107 | if ( id > 0 ) |
---|
| 108 | { |
---|
| 109 | // Append end; |
---|
| 110 | theSystem->SetParticipant ( new G4QMDParticipant ( (*it)->GetDefinition() , (*it)->Get4Momentum().v()/GeV , (*it)->GetPosition()/fermi ) ); |
---|
| 111 | et += (*it)->Get4Momentum().e()/GeV; |
---|
| 112 | if ( id > 1 ) |
---|
| 113 | { |
---|
| 114 | // std::cout << "NAGISA id >2; id= " << id << std::endl; |
---|
| 115 | } |
---|
| 116 | } |
---|
| 117 | id++; |
---|
| 118 | |
---|
| 119 | delete *it; |
---|
| 120 | } |
---|
| 121 | theMeanField->SetSystem ( theSystem ); |
---|
| 122 | i0 = id-1; // 0 enter to i |
---|
| 123 | } |
---|
| 124 | |
---|
| 125 | // EnergyCheck |
---|
| 126 | |
---|
| 127 | G4double epot = theMeanField->GetTotalPotential(); |
---|
| 128 | G4double efin = epot + et; |
---|
| 129 | //std::cout << std::abs ( eini - efin ) - epse << std::endl; |
---|
| 130 | // std::cout << std::abs ( eini - efin ) - epse*10 << std::endl; |
---|
| 131 | //071031 |
---|
| 132 | // *10 TK |
---|
| 133 | if ( std::abs ( eini - efin ) < epse*10 ) |
---|
| 134 | { |
---|
| 135 | // Energy OK |
---|
| 136 | // std::cout << "Decay Succeeded Energy OK" << std::endl; |
---|
| 137 | isThisEnergyOK = true; |
---|
| 138 | break; |
---|
| 139 | } |
---|
| 140 | else |
---|
| 141 | { |
---|
| 142 | for ( G4int i0i = 0 ; i0i < id-1 ; i0i++ ) |
---|
| 143 | { |
---|
| 144 | // std::cout << "Decay Energitically Blocked deleteing " << i0i+n0 << std::endl; |
---|
| 145 | theSystem->DeleteParticipant( i0i+n0 ); |
---|
| 146 | } |
---|
| 147 | } |
---|
| 148 | } |
---|
| 149 | |
---|
| 150 | // Pauli Check |
---|
| 151 | if ( isThisEnergyOK == true ) |
---|
| 152 | { |
---|
| 153 | // if ( theMeanField->IsPauliBlocked ( i ) != true ) |
---|
| 154 | { |
---|
| 155 | bool allOK = true; |
---|
| 156 | for ( G4int i0i = 0 ; i0i < i0 ; i0i++ ) |
---|
| 157 | { |
---|
| 158 | if ( theMeanField->IsPauliBlocked ( i0i+n0 ) == true ) |
---|
| 159 | { |
---|
| 160 | allOK = false; |
---|
| 161 | break; |
---|
| 162 | } |
---|
| 163 | } |
---|
| 164 | |
---|
| 165 | // if ( allOK ) std::cout << "Decay Succeeded" << std::endl; |
---|
| 166 | if ( allOK ) continue; //Do not Pauli Blocked |
---|
| 167 | } |
---|
| 168 | } |
---|
| 169 | // |
---|
| 170 | |
---|
| 171 | // std::cout << "Decay Blocked" << std::endl; |
---|
| 172 | theSystem->GetParticipant( i )->SetDefinition( pd0 ); |
---|
| 173 | theSystem->GetParticipant( i )->SetPosition( r0 ); |
---|
| 174 | theSystem->GetParticipant( i )->SetMomentum( p0 ); |
---|
| 175 | |
---|
| 176 | if ( isThisEnergyOK == true ) |
---|
| 177 | { |
---|
| 178 | for ( G4int i0i = 0 ; i0i < i0 ; i0i++ ) |
---|
| 179 | { |
---|
| 180 | // std::cout << "Decay Blocked deleteing " << i0i+n0 << std::endl; |
---|
| 181 | theSystem->DeleteParticipant( i0i+n0 ); |
---|
| 182 | } |
---|
| 183 | } |
---|
| 184 | |
---|
| 185 | } |
---|
| 186 | } |
---|
| 187 | //071101 |
---|
| 188 | |
---|
| 189 | |
---|
| 190 | n = theSystem->GetTotalNumberOfParticipant(); |
---|
| 191 | //std::cout << "Collision n " << n << std::endl; |
---|
| 192 | |
---|
| 193 | std::vector< G4bool > isCollided ( n , false ); |
---|
| 194 | |
---|
| 195 | for ( G4int i = 1 ; i < n ; i++ ) |
---|
| 196 | { |
---|
| 197 | |
---|
| 198 | //std::cout << "Collision i " << i << std::endl; |
---|
| 199 | |
---|
| 200 | G4ThreeVector ri = theSystem->GetParticipant( i )->GetPosition(); |
---|
| 201 | G4LorentzVector p4i = theSystem->GetParticipant( i )->Get4Momentum(); |
---|
| 202 | G4double rmi = theSystem->GetParticipant( i )->GetMass(); |
---|
| 203 | G4ParticleDefinition* pdi = theSystem->GetParticipant( i )->GetDefinition(); |
---|
| 204 | |
---|
| 205 | //std::cout << " p4i00 " << p4i << std::endl; |
---|
| 206 | for ( G4int j = 0 ; j < i ; j++ ) |
---|
| 207 | { |
---|
| 208 | // std::cout << "Collision " << i << " " << j << std::endl; |
---|
| 209 | |
---|
| 210 | /* |
---|
| 211 | std::cout << "Collision " << i << " " << theSystem->GetParticipant( i )->IsThisProjectile() << std::endl; |
---|
| 212 | std::cout << "Collision " << j << " " << theSystem->GetParticipant( j )->IsThisProjectile() << std::endl; |
---|
| 213 | std::cout << "Collision " << i << " " << theSystem->GetParticipant( i )->IsThisTarget() << std::endl; |
---|
| 214 | std::cout << "Collision " << j << " " << theSystem->GetParticipant( j )->IsThisTarget() << std::endl; |
---|
| 215 | */ |
---|
| 216 | |
---|
| 217 | // Only 1 Collision allowed for each particle in a time step. |
---|
| 218 | if ( isCollided[ i ] == true ) continue; |
---|
| 219 | if ( isCollided[ j ] == true ) continue; |
---|
| 220 | |
---|
| 221 | // Do not allow collision between nucleons in target/projectile til its first collision. |
---|
| 222 | if ( theSystem->GetParticipant( i )->IsThisProjectile() ) |
---|
| 223 | { |
---|
| 224 | if ( theSystem->GetParticipant( j )->IsThisProjectile() ) continue; |
---|
| 225 | } |
---|
| 226 | else if ( theSystem->GetParticipant( i )->IsThisTarget() ) |
---|
| 227 | { |
---|
| 228 | if ( theSystem->GetParticipant( j )->IsThisTarget() ) continue; |
---|
| 229 | } |
---|
| 230 | |
---|
| 231 | |
---|
| 232 | G4ThreeVector rj = theSystem->GetParticipant( j )->GetPosition(); |
---|
| 233 | G4LorentzVector p4j = theSystem->GetParticipant( j )->Get4Momentum(); |
---|
| 234 | G4double rmj = theSystem->GetParticipant( j )->GetMass(); |
---|
| 235 | G4ParticleDefinition* pdj = theSystem->GetParticipant( j )->GetDefinition(); |
---|
| 236 | |
---|
| 237 | G4double rr2 = theMeanField->GetRR2( i , j ); |
---|
| 238 | |
---|
| 239 | // Here we assume elab (beam momentum less than 5 GeV/n ) |
---|
| 240 | if ( rr2 > deltar*deltar ) continue; |
---|
| 241 | |
---|
| 242 | G4double s = (p4i+p4j)*(p4i+p4j); |
---|
| 243 | |
---|
| 244 | G4double srt = std::sqrt ( s ); |
---|
| 245 | |
---|
| 246 | G4double cutoff = 0.0; |
---|
| 247 | G4double bcmax = 0.0; |
---|
| 248 | G4double sig = 0.0; |
---|
| 249 | |
---|
| 250 | if ( rmi < 0.94 && rmj < 0.94 ) |
---|
| 251 | { |
---|
| 252 | // nucleon or pion case |
---|
| 253 | cutoff = rmi + rmj + 0.02; |
---|
| 254 | bcmax = bcmax0; |
---|
| 255 | sig = sig0; |
---|
| 256 | } |
---|
| 257 | else |
---|
| 258 | { |
---|
| 259 | cutoff = rmi + rmj; |
---|
| 260 | bcmax = bcmax1; |
---|
| 261 | sig = sig1; |
---|
| 262 | } |
---|
| 263 | |
---|
| 264 | //std::cout << "Collision cutoff " << i << " " << j << " " << cutoff << std::endl; |
---|
| 265 | if ( srt < cutoff ) continue; |
---|
| 266 | |
---|
| 267 | G4ThreeVector dr = ri - rj; |
---|
| 268 | G4double rsq = dr*dr; |
---|
| 269 | |
---|
| 270 | G4double pij = p4i*p4j; |
---|
| 271 | G4double pidr = p4i.vect()*dr; |
---|
| 272 | G4double pjdr = p4j.vect()*dr; |
---|
| 273 | |
---|
| 274 | G4double aij = 1.0 - ( rmi*rmj /pij ) * ( rmi*rmj /pij ); |
---|
| 275 | G4double bij = pidr / rmi - pjdr*rmi/pij; |
---|
| 276 | G4double cij = rsq + ( pidr / rmi ) * ( pidr / rmi ); |
---|
| 277 | G4double brel = std::sqrt ( std::abs ( cij - bij*bij/aij ) ); |
---|
| 278 | |
---|
| 279 | if ( brel > bcmax ) continue; |
---|
| 280 | //std::cout << "collisions3 " << std::endl; |
---|
| 281 | |
---|
| 282 | G4double bji = -pjdr/rmj + pidr * rmj /pij; |
---|
| 283 | |
---|
| 284 | G4double ti = ( pidr/rmi - bij / aij ) * p4i.e() / rmi; |
---|
| 285 | G4double tj = (-pjdr/rmj - bji / aij ) * p4j.e() / rmj; |
---|
| 286 | |
---|
| 287 | G4double deltaT = 0.0; |
---|
| 288 | deltaT = 1.0; // TK |
---|
| 289 | |
---|
| 290 | /* |
---|
| 291 | std::cout << "collisions4 p4i " << p4i << std::endl; |
---|
| 292 | std::cout << "collisions4 ri " << ri << std::endl; |
---|
| 293 | std::cout << "collisions4 p4j " << p4j << std::endl; |
---|
| 294 | std::cout << "collisions4 rj " << rj << std::endl; |
---|
| 295 | std::cout << "collisions4 dr " << dr << std::endl; |
---|
| 296 | std::cout << "collisions4 pij " << pij << std::endl; |
---|
| 297 | std::cout << "collisions4 aij " << aij << std::endl; |
---|
| 298 | std::cout << "collisions4 bij bji " << bij << " " << bji << std::endl; |
---|
| 299 | std::cout << "collisions4 pidr pjdr " << pidr << " " << pjdr << std::endl; |
---|
| 300 | std::cout << "collisions4 p4i.e() p4j.e() " << p4i.e() << " " << p4j.e() << std::endl; |
---|
| 301 | std::cout << "collisions4 rmi rmj " << rmi << " " << rmj << std::endl; |
---|
| 302 | std::cout << "collisions4 " << ti << " " << tj << std::endl; |
---|
| 303 | */ |
---|
| 304 | if ( std::abs ( ti + tj ) > deltaT ) continue; |
---|
| 305 | //std::cout << "collisions4 " << std::endl; |
---|
| 306 | |
---|
| 307 | G4ThreeVector beta = ( p4i + p4j ).boostVector(); |
---|
| 308 | |
---|
| 309 | G4LorentzVector p = p4i; |
---|
| 310 | G4LorentzVector p4icm = p.boost( p.findBoostToCM ( p4j ) ); |
---|
| 311 | G4ThreeVector pcm = p4icm.vect(); |
---|
| 312 | |
---|
| 313 | G4double prcm = pcm.mag(); |
---|
| 314 | |
---|
| 315 | if ( prcm <= 0.00001 ) continue; |
---|
| 316 | //std::cout << "collisions5 " << std::endl; |
---|
| 317 | |
---|
| 318 | G4bool energetically_forbidden = !( CalFinalStateOfTheBinaryCollision ( i , j ) ); // Use Geant4 Collision Library |
---|
| 319 | //G4bool energetically_forbidden = !( CalFinalStateOfTheBinaryCollisionJQMD ( sig , cutoff , pcm , prcm , srt, beta , gamma , i , j ) ); // JQMD Elastic |
---|
| 320 | |
---|
| 321 | G4bool pauli_blocked = false; |
---|
| 322 | if ( energetically_forbidden != true ) |
---|
| 323 | { |
---|
| 324 | if ( theMeanField->IsPauliBlocked ( i ) == true || theMeanField->IsPauliBlocked ( j ) == true ) |
---|
| 325 | { |
---|
| 326 | pauli_blocked = true; |
---|
| 327 | //std::cout << "G4QMDRESULT Collsion Pauli Blocked " << std::endl; |
---|
| 328 | } |
---|
| 329 | } |
---|
| 330 | else |
---|
| 331 | { |
---|
| 332 | //std::cout << "G4QMDRESULT Collsion Blocked " << std::endl; |
---|
| 333 | } |
---|
| 334 | |
---|
| 335 | /* |
---|
| 336 | std::cout << "G4QMDRESULT Collsion initial p4 i and j " |
---|
| 337 | << p4i << " " << p4j |
---|
| 338 | << std::endl; |
---|
| 339 | */ |
---|
| 340 | |
---|
| 341 | |
---|
| 342 | if ( energetically_forbidden == true || pauli_blocked == true ) |
---|
| 343 | { |
---|
| 344 | // Collsion not allowed then re enter orginal participants |
---|
| 345 | // Now only momentum, becasuse we only consider elastic scattering of nucleons |
---|
| 346 | |
---|
| 347 | theSystem->GetParticipant( i )->SetMomentum( p4i.vect() ); |
---|
| 348 | theSystem->GetParticipant( i )->SetDefinition( pdi ); |
---|
| 349 | theSystem->GetParticipant( i )->SetPosition( ri ); |
---|
| 350 | theSystem->GetParticipant( j )->SetMomentum( p4j.vect() ); |
---|
| 351 | theSystem->GetParticipant( j )->SetDefinition( pdj ); |
---|
| 352 | theSystem->GetParticipant( j )->SetPosition( rj ); |
---|
| 353 | |
---|
| 354 | } |
---|
| 355 | else |
---|
| 356 | { |
---|
| 357 | // Collsion allowed (really happened) |
---|
| 358 | |
---|
| 359 | // Unset Projectile/Target flag |
---|
| 360 | theSystem->GetParticipant( i )->UnsetInitialMark(); |
---|
| 361 | theSystem->GetParticipant( j )->UnsetInitialMark(); |
---|
| 362 | |
---|
| 363 | isCollided[ i ] = true; |
---|
| 364 | isCollided[ j ] = true; |
---|
| 365 | |
---|
| 366 | theSystem->IncrementCollisionCounter(); |
---|
| 367 | |
---|
| 368 | /* |
---|
| 369 | std::cout << "G4QMDRESULT Collsion Really Happened between " |
---|
| 370 | << i << " and " << j |
---|
| 371 | << std::endl; |
---|
| 372 | std::cout << "G4QMDRESULT Collsion initial p4 i and j " |
---|
| 373 | << p4i << " " << p4j |
---|
| 374 | << std::endl; |
---|
| 375 | std::cout << "G4QMDRESULT Collsion after p4 i and j " |
---|
| 376 | << theSystem->GetParticipant( i )->Get4Momentum() |
---|
| 377 | << " " |
---|
| 378 | << theSystem->GetParticipant( j )->Get4Momentum() |
---|
| 379 | << std::endl; |
---|
| 380 | std::cout << "G4QMDRESULT Collsion Diff " |
---|
| 381 | << p4i + p4j - theSystem->GetParticipant( i )->Get4Momentum() - theSystem->GetParticipant( j )->Get4Momentum() |
---|
| 382 | << std::endl; |
---|
| 383 | std::cout << "G4QMDRESULT Collsion initial r i and j " |
---|
| 384 | << ri << " " << rj |
---|
| 385 | << std::endl; |
---|
| 386 | std::cout << "G4QMDRESULT Collsion after r i and j " |
---|
| 387 | << theSystem->GetParticipant( i )->GetPosition() |
---|
| 388 | << " " |
---|
| 389 | << theSystem->GetParticipant( j )->GetPosition() |
---|
| 390 | << std::endl; |
---|
| 391 | */ |
---|
| 392 | |
---|
| 393 | } |
---|
| 394 | |
---|
| 395 | // theMeanField |
---|
| 396 | |
---|
| 397 | } |
---|
| 398 | } |
---|
| 399 | |
---|
| 400 | //071106 |
---|
| 401 | n = theSystem->GetTotalNumberOfParticipant(); |
---|
| 402 | G4bool isThisModefied = false; |
---|
| 403 | for ( G4int i = 0 ; i < n ; i++ ) |
---|
| 404 | { |
---|
| 405 | if ( theSystem->GetParticipant( i )->GetDefinition() == G4PionZero::PionZero() ) |
---|
| 406 | { |
---|
| 407 | if ( theSystem->GetParticipant( i )->GetPosition().mag() > 1.0e9 ) |
---|
| 408 | { |
---|
| 409 | // std::cout << "Deleting " << i << " " << theSystem->GetParticipant( i )->GetPosition().mag() << std::endl; |
---|
| 410 | theSystem->DeleteParticipant( i ); |
---|
| 411 | isThisModefied = true; |
---|
| 412 | } |
---|
| 413 | } |
---|
| 414 | } |
---|
| 415 | if ( isThisModefied == true ) theMeanField->SetSystem ( theSystem ); |
---|
| 416 | //071106 |
---|
| 417 | |
---|
| 418 | } |
---|
| 419 | |
---|
| 420 | |
---|
| 421 | |
---|
| 422 | G4bool G4QMDCollision::CalFinalStateOfTheBinaryCollision( G4int i , G4int j ) |
---|
| 423 | { |
---|
| 424 | |
---|
| 425 | //G4cout << "CalFinalStateOfTheBinaryCollision " << G4endl; |
---|
| 426 | |
---|
| 427 | G4bool result = true; |
---|
| 428 | |
---|
| 429 | G4LorentzVector p4i = theSystem->GetParticipant( i )->Get4Momentum(); |
---|
| 430 | G4LorentzVector p4j = theSystem->GetParticipant( j )->Get4Momentum(); |
---|
| 431 | |
---|
| 432 | //071031 |
---|
| 433 | // will use KineticTrack |
---|
| 434 | G4LorentzVector p4ix = p4i*GeV; |
---|
| 435 | G4LorentzVector p4jx = p4j*GeV; |
---|
| 436 | G4ThreeVector rix = (theSystem->GetParticipant( i )->GetPosition())*fermi; |
---|
| 437 | G4ThreeVector rjx = (theSystem->GetParticipant( j )->GetPosition())*fermi; |
---|
| 438 | //071031 |
---|
| 439 | |
---|
| 440 | G4double epot = theMeanField->GetTotalPotential(); |
---|
| 441 | |
---|
| 442 | G4double eini = epot + p4i.e() + p4j.e(); |
---|
| 443 | |
---|
| 444 | |
---|
| 445 | //071031 |
---|
| 446 | G4ParticleDefinition* pdi0 =theSystem->GetParticipant( i )->GetDefinition(); |
---|
| 447 | G4ParticleDefinition* pdj0 =theSystem->GetParticipant( j )->GetDefinition(); |
---|
| 448 | G4ThreeVector ri0 =(theSystem->GetParticipant( i )->GetPosition())*fermi; |
---|
| 449 | G4ThreeVector rj0 =(theSystem->GetParticipant( j )->GetPosition())*fermi; |
---|
| 450 | |
---|
| 451 | for ( G4int iitry = 0 ; iitry < 4 ; iitry++ ) |
---|
| 452 | { |
---|
| 453 | |
---|
| 454 | G4KineticTrack kt1( pdi0 , 0.0 , ri0 , p4ix ); |
---|
| 455 | G4KineticTrack kt2( pdj0 , 0.0 , rj0 , p4jx ); |
---|
| 456 | G4LorentzVector p4ix_new; |
---|
| 457 | G4LorentzVector p4jx_new; |
---|
| 458 | G4KineticTrackVector* secs = NULL; |
---|
| 459 | secs = theScatterer->Scatter( kt1 , kt2 ); |
---|
| 460 | |
---|
| 461 | //std::cout << "G4QMDSCATTERER BEFORE " << kt1.GetDefinition()->GetParticleName() << " " << kt1.Get4Momentum()/GeV << " " << kt1.GetPosition()/fermi << std::endl; |
---|
| 462 | //std::cout << "G4QMDSCATTERER BEFORE " << kt2.GetDefinition()->GetParticleName() << " " << kt2.Get4Momentum()/GeV << " " << kt2.GetPosition()/fermi << std::endl; |
---|
| 463 | //std::cout << "THESCATTERER " << theScatterer->GetCrossSection ( kt1 , kt2 )/millibarn << " " << elastic << " " << sig << std::endl; |
---|
| 464 | if ( secs ) |
---|
| 465 | { |
---|
| 466 | G4int iti = 0; |
---|
| 467 | if ( secs->size() == 2 ) |
---|
| 468 | { |
---|
| 469 | for ( G4KineticTrackVector::iterator it |
---|
| 470 | = secs->begin() ; it != secs->end() ; it++ ) |
---|
| 471 | { |
---|
| 472 | if ( iti == 0 ) |
---|
| 473 | { |
---|
| 474 | theSystem->GetParticipant( i )->SetDefinition( (*it)->GetDefinition() ); |
---|
| 475 | p4ix_new = (*it)->Get4Momentum()/GeV; |
---|
| 476 | //std::cout << "THESCATTERER " << (*it)->GetDefinition()->GetParticleName() << std::endl; |
---|
| 477 | theSystem->GetParticipant( i )->SetMomentum( p4ix_new.v() ); |
---|
| 478 | } |
---|
| 479 | if ( iti == 1 ) |
---|
| 480 | { |
---|
| 481 | theSystem->GetParticipant( j )->SetDefinition( (*it)->GetDefinition() ); |
---|
| 482 | p4jx_new = (*it)->Get4Momentum()/GeV; |
---|
| 483 | //std::cout << "THESCATTERER " << p4jx_new.e()-p4jx_new.m() << std::endl; |
---|
| 484 | theSystem->GetParticipant( j )->SetMomentum( p4jx_new.v() ); |
---|
| 485 | } |
---|
| 486 | //std::cout << "G4QMDSCATTERER AFTER " << (*it)->GetDefinition()->GetParticleName() << " " << (*it)->Get4Momentum()/GeV << std::endl; |
---|
| 487 | iti++; |
---|
| 488 | } |
---|
| 489 | } |
---|
| 490 | else |
---|
| 491 | { |
---|
| 492 | //std::cout << "NAGISA pion absrorption " << secs->front()->GetDefinition()->GetParticleName() << std::endl; |
---|
| 493 | //secs->front()->Decay(); |
---|
| 494 | theSystem->GetParticipant( i )->SetDefinition( secs->front()->GetDefinition() ); |
---|
| 495 | p4ix_new = secs->front()->Get4Momentum()/GeV; |
---|
| 496 | theSystem->GetParticipant( i )->SetMomentum( p4ix_new.v() ); |
---|
| 497 | |
---|
| 498 | //std::cout << "THESCATTERER " << (*it)->GetDefinition()->GetParticleName() << std::endl; |
---|
| 499 | p4jx_new( 0 ); |
---|
| 500 | //theSystem->GetParticipant( j )->SetDefinition( G4Gamma::Gamma() ); |
---|
| 501 | //theSystem->GetParticipant( j )->SetDefinition( G4Neutron::Neutron() ); |
---|
| 502 | theSystem->GetParticipant( j )->SetDefinition( G4PionZero::PionZero() ); |
---|
| 503 | theSystem->GetParticipant( j )->SetMomentum( G4ThreeVector( G4UniformRand() )*eV ); |
---|
| 504 | theSystem->GetParticipant( j )->SetPosition( G4ThreeVector( 1000, 1000, 1000 )*km ); |
---|
| 505 | |
---|
| 506 | } |
---|
| 507 | |
---|
| 508 | if ( secs->size() > 2 ) std::cout << "NAGISA secs size > 2; " << secs->size() << std::endl; |
---|
| 509 | |
---|
| 510 | // deleteing KineticTrack |
---|
| 511 | for ( G4KineticTrackVector::iterator it |
---|
| 512 | = secs->begin() ; it != secs->end() ; it++ ) |
---|
| 513 | { |
---|
| 514 | delete *it; |
---|
| 515 | } |
---|
| 516 | } |
---|
| 517 | //071031 |
---|
| 518 | |
---|
| 519 | theMeanField->Cal2BodyQuantities( i ); |
---|
| 520 | theMeanField->Cal2BodyQuantities( j ); |
---|
| 521 | |
---|
| 522 | epot = theMeanField->GetTotalPotential(); |
---|
| 523 | |
---|
| 524 | G4double efin = epot + p4ix_new.e() + p4jx_new.e(); |
---|
| 525 | |
---|
| 526 | //std::cout << "Collision NEW epot " << i << " " << j << " " << epot << " " << std::abs ( eini - efin ) - epse << std::endl; |
---|
| 527 | |
---|
| 528 | /* |
---|
| 529 | std::cout << "Collision efin " << i << " " << j << " " << efin << std::endl; |
---|
| 530 | std::cout << "Collision " << i << " " << j << " " << std::abs ( eini - efin ) << " " << epse << std::endl; |
---|
| 531 | std::cout << "Collision " << std::abs ( eini - efin ) << " " << epse << std::endl; |
---|
| 532 | */ |
---|
| 533 | |
---|
| 534 | //071031 |
---|
| 535 | if ( std::abs ( eini - efin ) < epse ) |
---|
| 536 | { |
---|
| 537 | // Collison OK |
---|
| 538 | //std::cout << "collisions6" << std::endl; |
---|
| 539 | //std::cout << "collisions before " << p4i << " " << p4j << std::endl; |
---|
| 540 | //std::cout << "collisions after " << theSystem->GetParticipant( i )->Get4Momentum() << " " << theSystem->GetParticipant( j )->Get4Momentum() << std::endl; |
---|
| 541 | //std::cout << "collisions dif " << ( p4i + p4j ) - ( theSystem->GetParticipant( i )->Get4Momentum() + theSystem->GetParticipant( j )->Get4Momentum() ) << std::endl; |
---|
| 542 | //std::cout << "collisions before " << rix/fermi << " " << rjx/fermi << std::endl; |
---|
| 543 | //std::cout << "collisions after " << theSystem->GetParticipant( i )->GetPosition() << " " << theSystem->GetParticipant( j )->GetPosition() << std::endl; |
---|
| 544 | } |
---|
| 545 | //071031 |
---|
| 546 | |
---|
| 547 | if ( std::abs ( eini - efin ) < epse ) return result; // Collison OK |
---|
| 548 | |
---|
| 549 | } |
---|
| 550 | |
---|
| 551 | // Energetically forbidden collision |
---|
| 552 | result = false; |
---|
| 553 | |
---|
| 554 | return result; |
---|
| 555 | |
---|
| 556 | } |
---|
| 557 | |
---|
| 558 | |
---|
| 559 | |
---|
| 560 | G4bool G4QMDCollision::CalFinalStateOfTheBinaryCollisionJQMD( G4double sig , G4double cutoff , G4ThreeVector pcm , G4double prcm , G4double srt , G4ThreeVector beta , G4double gamma , G4int i , G4int j ) |
---|
| 561 | { |
---|
| 562 | |
---|
| 563 | //G4cout << "CalFinalStateOfTheBinaryCollisionJQMD" << G4endl; |
---|
| 564 | |
---|
| 565 | G4bool result = true; |
---|
| 566 | |
---|
| 567 | G4LorentzVector p4i = theSystem->GetParticipant( i )->Get4Momentum(); |
---|
| 568 | G4double rmi = theSystem->GetParticipant( i )->GetMass(); |
---|
| 569 | G4int zi = theSystem->GetParticipant( i )->GetChargeInUnitOfEplus(); |
---|
| 570 | |
---|
| 571 | G4LorentzVector p4j = theSystem->GetParticipant( j )->Get4Momentum(); |
---|
| 572 | G4double rmj = theSystem->GetParticipant( j )->GetMass(); |
---|
| 573 | G4int zj = theSystem->GetParticipant( j )->GetChargeInUnitOfEplus(); |
---|
| 574 | |
---|
| 575 | G4double pr = prcm; |
---|
| 576 | |
---|
| 577 | G4double c2 = pcm.z()/pr; |
---|
| 578 | |
---|
| 579 | G4double csrt = srt - cutoff; |
---|
| 580 | |
---|
| 581 | //G4double pri = prcm; |
---|
| 582 | //G4double prf = sqrt ( 0.25 * srt*srt -rm2 ); |
---|
| 583 | |
---|
| 584 | G4double asrt = srt - rmi - rmj; |
---|
| 585 | G4double pra = prcm; |
---|
| 586 | |
---|
| 587 | |
---|
| 588 | |
---|
| 589 | G4double elastic = 0.0; |
---|
| 590 | |
---|
| 591 | if ( zi == zj ) |
---|
| 592 | { |
---|
| 593 | if ( csrt < 0.4286 ) |
---|
| 594 | { |
---|
| 595 | elastic = 35.0 / ( 1. + csrt * 100.0 ) + 20.0; |
---|
| 596 | } |
---|
| 597 | else |
---|
| 598 | { |
---|
| 599 | elastic = ( - std::atan( ( csrt - 0.4286 ) * 1.5 - 0.8 ) |
---|
| 600 | * 2. / pi + 1.0 ) * 9.65 + 7.0; |
---|
| 601 | } |
---|
| 602 | } |
---|
| 603 | else |
---|
| 604 | { |
---|
| 605 | if ( csrt < 0.4286 ) |
---|
| 606 | { |
---|
| 607 | elastic = 28.0 / ( 1. + csrt * 100.0 ) + 27.0; |
---|
| 608 | } |
---|
| 609 | else |
---|
| 610 | { |
---|
| 611 | elastic = ( - std::atan( ( csrt - 0.4286 ) * 1.5 - 0.8 ) |
---|
| 612 | * 2. / pi + 1.0 ) * 12.34 + 10.0; |
---|
| 613 | } |
---|
| 614 | } |
---|
| 615 | |
---|
| 616 | // std::cout << "Collision csrt " << i << " " << j << " " << csrt << std::endl; |
---|
| 617 | // std::cout << "Collision elstic " << i << " " << j << " " << elastic << std::endl; |
---|
| 618 | |
---|
| 619 | |
---|
| 620 | // std::cout << "Collision sig " << i << " " << j << " " << sig << std::endl; |
---|
| 621 | if ( G4UniformRand() > elastic / sig ) |
---|
| 622 | { |
---|
| 623 | //std::cout << "Inelastic " << std::endl; |
---|
| 624 | //std::cout << "elastic/sig " << elastic/sig << std::endl; |
---|
| 625 | return result; |
---|
| 626 | } |
---|
| 627 | else |
---|
| 628 | { |
---|
| 629 | //std::cout << "Elastic " << std::endl; |
---|
| 630 | } |
---|
| 631 | // std::cout << "Collision ELSTIC " << i << " " << j << std::endl; |
---|
| 632 | |
---|
| 633 | |
---|
| 634 | G4double as = std::pow ( 3.65 * asrt , 6 ); |
---|
| 635 | G4double a = 6.0 * as / (1.0 + as); |
---|
| 636 | G4double ta = -2.0 * pra*pra; |
---|
| 637 | G4double x = G4UniformRand(); |
---|
| 638 | G4double t1 = std::log( (1-x) * std::exp(2.*a*ta) + x ) / a; |
---|
| 639 | G4double c1 = 1.0 - t1/ta; |
---|
| 640 | |
---|
| 641 | if( std::abs(c1) > 1.0 ) c1 = 2.0 * x - 1.0; |
---|
| 642 | |
---|
| 643 | /* |
---|
| 644 | std::cout << "Collision as " << i << " " << j << " " << as << std::endl; |
---|
| 645 | std::cout << "Collision a " << i << " " << j << " " << a << std::endl; |
---|
| 646 | std::cout << "Collision ta " << i << " " << j << " " << ta << std::endl; |
---|
| 647 | std::cout << "Collision x " << i << " " << j << " " << x << std::endl; |
---|
| 648 | std::cout << "Collision t1 " << i << " " << j << " " << t1 << std::endl; |
---|
| 649 | std::cout << "Collision c1 " << i << " " << j << " " << c1 << std::endl; |
---|
| 650 | */ |
---|
| 651 | t1 = 2.0*pi*G4UniformRand(); |
---|
| 652 | // std::cout << "Collision t1 " << i << " " << j << " " << t1 << std::endl; |
---|
| 653 | G4double t2 = 0.0; |
---|
| 654 | if ( pcm.x() == 0.0 && pcm.y() == 0 ) |
---|
| 655 | { |
---|
| 656 | t2 = 0.0; |
---|
| 657 | } |
---|
| 658 | else |
---|
| 659 | { |
---|
| 660 | t2 = std::atan2( pcm.y() , pcm.x() ); |
---|
| 661 | } |
---|
| 662 | // std::cout << "Collision t2 " << i << " " << j << " " << t2 << std::endl; |
---|
| 663 | |
---|
| 664 | G4double s1 = std::sqrt ( 1.0 - c1*c1 ); |
---|
| 665 | G4double s2 = std::sqrt ( 1.0 - c2*c2 ); |
---|
| 666 | |
---|
| 667 | G4double ct1 = std::cos(t1); |
---|
| 668 | G4double st1 = std::sin(t1); |
---|
| 669 | |
---|
| 670 | G4double ct2 = std::cos(t2); |
---|
| 671 | G4double st2 = std::sin(t2); |
---|
| 672 | |
---|
| 673 | G4double ss = c2*s1*ct1 + s2*c1; |
---|
| 674 | |
---|
| 675 | pcm.setX( pr * ( ss*ct2 - s1*st1*st2) ); |
---|
| 676 | pcm.setY( pr * ( ss*st2 + s1*st1*ct2) ); |
---|
| 677 | pcm.setZ( pr * ( c1*c2 - s1*s2*ct1) ); |
---|
| 678 | |
---|
| 679 | // std::cout << "Collision pcm " << i << " " << j << " " << pcm << std::endl; |
---|
| 680 | |
---|
| 681 | G4double epot = theMeanField->GetTotalPotential(); |
---|
| 682 | |
---|
| 683 | G4double eini = epot + p4i.e() + p4j.e(); |
---|
| 684 | G4double etwo = p4i.e() + p4j.e(); |
---|
| 685 | |
---|
| 686 | /* |
---|
| 687 | std::cout << "Collision epot " << i << " " << j << " " << epot << std::endl; |
---|
| 688 | std::cout << "Collision eini " << i << " " << j << " " << eini << std::endl; |
---|
| 689 | std::cout << "Collision etwo " << i << " " << j << " " << etwo << std::endl; |
---|
| 690 | */ |
---|
| 691 | |
---|
| 692 | |
---|
| 693 | for ( G4int itry = 0 ; itry < 4 ; itry++ ) |
---|
| 694 | { |
---|
| 695 | |
---|
| 696 | G4double eicm = std::sqrt ( rmi*rmi + pcm*pcm ); |
---|
| 697 | G4double pibeta = pcm*beta; |
---|
| 698 | |
---|
| 699 | G4double trans = gamma * ( gamma * pibeta / ( gamma + 1 ) + eicm ); |
---|
| 700 | |
---|
| 701 | G4ThreeVector pi_new = beta*trans + pcm; |
---|
| 702 | |
---|
| 703 | G4double ejcm = std::sqrt ( rmj*rmj + pcm*pcm ); |
---|
| 704 | trans = gamma * ( gamma * pibeta / ( gamma + 1 ) + ejcm ); |
---|
| 705 | |
---|
| 706 | G4ThreeVector pj_new = beta*trans - pcm; |
---|
| 707 | |
---|
| 708 | // |
---|
| 709 | // Delete old |
---|
| 710 | // Add new Particitipants |
---|
| 711 | // |
---|
| 712 | // Now only change momentum ( Beacuse we only have elastic sctter of nucleon |
---|
| 713 | // In future Definition also will be change |
---|
| 714 | // |
---|
| 715 | |
---|
| 716 | theSystem->GetParticipant( i )->SetMomentum( pi_new ); |
---|
| 717 | theSystem->GetParticipant( j )->SetMomentum( pj_new ); |
---|
| 718 | |
---|
| 719 | G4double pi_new_e = (theSystem->GetParticipant( i )->Get4Momentum()).e(); |
---|
| 720 | G4double pj_new_e = (theSystem->GetParticipant( j )->Get4Momentum()).e(); |
---|
| 721 | |
---|
| 722 | theMeanField->Cal2BodyQuantities( i ); |
---|
| 723 | theMeanField->Cal2BodyQuantities( j ); |
---|
| 724 | |
---|
| 725 | epot = theMeanField->GetTotalPotential(); |
---|
| 726 | |
---|
| 727 | G4double efin = epot + pi_new_e + pj_new_e ; |
---|
| 728 | |
---|
| 729 | //std::cout << "Collision NEW epot " << i << " " << j << " " << epot << " " << std::abs ( eini - efin ) - epse << std::endl; |
---|
| 730 | /* |
---|
| 731 | std::cout << "Collision efin " << i << " " << j << " " << efin << std::endl; |
---|
| 732 | std::cout << "Collision " << i << " " << j << " " << std::abs ( eini - efin ) << " " << epse << std::endl; |
---|
| 733 | std::cout << "Collision " << std::abs ( eini - efin ) << " " << epse << std::endl; |
---|
| 734 | */ |
---|
| 735 | |
---|
| 736 | //071031 |
---|
| 737 | if ( std::abs ( eini - efin ) < epse ) |
---|
| 738 | { |
---|
| 739 | // Collison OK |
---|
| 740 | //std::cout << "collisions6" << std::endl; |
---|
| 741 | //std::cout << "collisions before " << p4i << " " << p4j << std::endl; |
---|
| 742 | //std::cout << "collisions after " << theSystem->GetParticipant( i )->Get4Momentum() << " " << theSystem->GetParticipant( j )->Get4Momentum() << std::endl; |
---|
| 743 | //std::cout << "collisions dif " << ( p4i + p4j ) - ( theSystem->GetParticipant( i )->Get4Momentum() + theSystem->GetParticipant( j )->Get4Momentum() ) << std::endl; |
---|
| 744 | //std::cout << "collisions before " << rix/fermi << " " << rjx/fermi << std::endl; |
---|
| 745 | //std::cout << "collisions after " << theSystem->GetParticipant( i )->GetPosition() << " " << theSystem->GetParticipant( j )->GetPosition() << std::endl; |
---|
| 746 | } |
---|
| 747 | //071031 |
---|
| 748 | |
---|
| 749 | if ( std::abs ( eini - efin ) < epse ) return result; // Collison OK |
---|
| 750 | |
---|
| 751 | G4double cona = ( eini - efin + etwo ) / gamma; |
---|
| 752 | G4double fac2 = 1.0 / ( 4.0 * cona*cona * pr*pr ) * |
---|
| 753 | ( ( cona*cona - ( rmi*rmi + rmj*rmj ) )*( cona*cona - ( rmi*rmi + rmj*rmj ) ) |
---|
| 754 | - 4.0 * rmi*rmi * rmj*rmj ); |
---|
| 755 | |
---|
| 756 | if ( fac2 > 0 ) |
---|
| 757 | { |
---|
| 758 | G4double fact = std::sqrt ( fac2 ); |
---|
| 759 | pcm = fact*pcm; |
---|
| 760 | } |
---|
| 761 | |
---|
| 762 | |
---|
| 763 | } |
---|
| 764 | |
---|
| 765 | // Energetically forbidden collision |
---|
| 766 | result = false; |
---|
| 767 | |
---|
| 768 | return result; |
---|
| 769 | |
---|
| 770 | } |
---|