[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 | // |
---|
[1055] | 27 | // Created: H.-P. Wellicsh: GHAD model wrapper for the CHIPS model (1997) |
---|
[819] | 28 | // 16.01.08 V.Ivanchenko use initialization similar to other CHIPS models |
---|
| 29 | // |
---|
| 30 | |
---|
| 31 | #include "G4ChiralInvariantPhaseSpace.hh" |
---|
| 32 | #include "G4ParticleTable.hh" |
---|
| 33 | #include "G4LorentzVector.hh" |
---|
| 34 | #include "G4DynamicParticle.hh" |
---|
| 35 | #include "G4IonTable.hh" |
---|
| 36 | #include "G4Neutron.hh" |
---|
| 37 | |
---|
| 38 | G4ChiralInvariantPhaseSpace::G4ChiralInvariantPhaseSpace() |
---|
| 39 | {} |
---|
| 40 | |
---|
| 41 | G4ChiralInvariantPhaseSpace::~G4ChiralInvariantPhaseSpace() |
---|
| 42 | {} |
---|
| 43 | |
---|
| 44 | G4HadFinalState * G4ChiralInvariantPhaseSpace::ApplyYourself( |
---|
| 45 | const G4HadProjectile& aTrack, |
---|
| 46 | G4Nucleus& aTargetNucleus, |
---|
| 47 | G4HadFinalState * aChange) |
---|
| 48 | { |
---|
| 49 | G4HadFinalState * aResult; |
---|
| 50 | if(aChange != 0) |
---|
| 51 | { |
---|
| 52 | aResult = aChange; |
---|
| 53 | } |
---|
| 54 | else |
---|
| 55 | { |
---|
| 56 | aResult = & theResult; |
---|
| 57 | aResult->Clear(); |
---|
| 58 | aResult->SetStatusChange(stopAndKill); |
---|
| 59 | } |
---|
| 60 | //projectile properties needed in constructor of quasmon |
---|
| 61 | G4LorentzVector proj4Mom; |
---|
| 62 | proj4Mom = aTrack.Get4Momentum(); |
---|
| 63 | G4int projectilePDGCode = aTrack.GetDefinition() |
---|
| 64 | ->GetPDGEncoding(); |
---|
| 65 | |
---|
| 66 | //target properties needed in constructor of quasmon |
---|
| 67 | G4int targetZ = G4int(aTargetNucleus.GetZ()+0.5); |
---|
| 68 | G4int targetA = G4int(aTargetNucleus.GetN()+0.5); |
---|
| 69 | G4int targetPDGCode = 90000000 + 1000*targetZ + (targetA-targetZ); |
---|
| 70 | // NOT NECESSARY ______________ |
---|
| 71 | G4double targetMass = G4ParticleTable::GetParticleTable()->GetIonTable() |
---|
| 72 | ->GetIonMass(targetZ, targetA); |
---|
| 73 | G4LorentzVector targ4Mom(0.,0.,0.,targetMass); |
---|
| 74 | // END OF NOT NECESSARY^^^^^^^^ |
---|
| 75 | |
---|
| 76 | // V.Ivanchenko set the same value as for other CHIPS models |
---|
| 77 | G4int nop = 152; |
---|
| 78 | //G4int nop = 164; // nuclear clusters up to A=21 |
---|
| 79 | G4double fractionOfSingleQuasiFreeNucleons = 0.4; |
---|
| 80 | G4double fractionOfPairedQuasiFreeNucleons = 0.0; |
---|
| 81 | if(targetA>27) fractionOfPairedQuasiFreeNucleons = 0.04; |
---|
| 82 | G4double clusteringCoefficient = 4.; |
---|
| 83 | G4double temperature = 180.; |
---|
| 84 | G4double halfTheStrangenessOfSee = 0.1; // = s/d = s/u |
---|
| 85 | G4double etaToEtaPrime = 0.3; |
---|
| 86 | |
---|
| 87 | // construct and fragment the quasmon |
---|
| 88 | //G4QCHIPSWorld aWorld(nop); // Create CHIPS World of nop particles |
---|
| 89 | G4QCHIPSWorld::Get()->GetParticles(nop); // Create CHIPS World of nop particles |
---|
| 90 | G4QNucleus::SetParameters(fractionOfSingleQuasiFreeNucleons, |
---|
| 91 | fractionOfPairedQuasiFreeNucleons, |
---|
| 92 | clusteringCoefficient); |
---|
| 93 | G4Quasmon::SetParameters(temperature, |
---|
| 94 | halfTheStrangenessOfSee, |
---|
| 95 | etaToEtaPrime); |
---|
| 96 | // G4QEnvironment::SetParameters(solidAngle); |
---|
[1055] | 97 | // G4cout << "Input info "<< projectilePDGCode << " " |
---|
| 98 | // << targetPDGCode <<" " |
---|
| 99 | // << 1./MeV*proj4Mom<<" " |
---|
| 100 | // << 1./MeV*targ4Mom << " " |
---|
| 101 | // << nop << G4endl; |
---|
[819] | 102 | G4QHadronVector projHV; |
---|
| 103 | G4QHadron* iH = new G4QHadron(projectilePDGCode, 1./MeV*proj4Mom); |
---|
| 104 | projHV.push_back(iH); |
---|
| 105 | G4QEnvironment* pan= new G4QEnvironment(projHV, targetPDGCode); |
---|
[1055] | 106 | //G4Quasmon* pan= new G4Quasmon(projectilePDGCode, targetPDGCode, |
---|
| 107 | // 1./MeV*proj4Mom, 1./MeV*targ4Mom, nop); |
---|
[819] | 108 | G4QHadronVector* output=0; |
---|
| 109 | try |
---|
| 110 | { |
---|
| 111 | output = pan->Fragment(); |
---|
| 112 | } |
---|
| 113 | catch(G4HadronicException & aR) |
---|
| 114 | { |
---|
| 115 | G4cerr << "Exception thrown passing through G4ChiralInvariantPhaseSpace "<<G4endl; |
---|
| 116 | G4cerr << " targetPDGCode = "<< targetPDGCode <<G4endl; |
---|
| 117 | G4cerr << " Dumping the information in the pojectile list"<<G4endl; |
---|
| 118 | for(size_t i=0; i< projHV.size(); i++) |
---|
| 119 | { |
---|
| 120 | G4cerr <<" Incoming 4-momentum and PDG code of "<<i<<"'th hadron: " |
---|
| 121 | <<" "<< projHV[i]->Get4Momentum()<<" "<<projHV[i]->GetPDGCode()<<G4endl; |
---|
| 122 | } |
---|
| 123 | throw; |
---|
| 124 | } |
---|
| 125 | std::for_each(projHV.begin(), projHV.end(), DeleteQHadron()); |
---|
| 126 | projHV.clear(); |
---|
| 127 | delete pan; |
---|
| 128 | |
---|
| 129 | // Fill the particle change. |
---|
| 130 | G4DynamicParticle * theSec; |
---|
| 131 | #ifdef CHIPSdebug |
---|
| 132 | G4cout << "G4ChiralInvariantPhaseSpace: NEW EVENT #ofHadrons=" |
---|
| 133 | <<output->size()<<G4endl; |
---|
| 134 | #endif |
---|
| 135 | unsigned int particle; |
---|
| 136 | for( particle = 0; particle < output->size(); particle++) |
---|
| 137 | { |
---|
| 138 | if(output->operator[](particle)->GetNFragments() != 0) |
---|
| 139 | { |
---|
| 140 | delete output->operator[](particle); |
---|
| 141 | continue; |
---|
| 142 | } |
---|
| 143 | theSec = new G4DynamicParticle; |
---|
| 144 | G4int pdgCode = output->operator[](particle)->GetPDGCode(); |
---|
| 145 | #ifdef CHIPSdebug |
---|
| 146 | G4cout << "G4ChiralInvariantPhaseSpace: h#"<<particle |
---|
| 147 | <<", PDG="<<pdgCode<<G4endl; |
---|
| 148 | #endif |
---|
| 149 | G4ParticleDefinition * theDefinition; |
---|
| 150 | // Note that I still have to take care of strange nuclei |
---|
| 151 | // For this I need the mass calculation, and a changed interface |
---|
| 152 | // for ion-tablel ==> work for Hisaya @@@@@@@ |
---|
| 153 | // Then I can sort out the pdgCode. I also need a decau process |
---|
| 154 | // for strange nuclei; may be another chips interface |
---|
| 155 | if(pdgCode>90000000) |
---|
| 156 | { |
---|
| 157 | G4int aZ = (pdgCode-90000000)/1000; |
---|
| 158 | G4int anN = pdgCode-90000000-1000*aZ; |
---|
| 159 | theDefinition = G4ParticleTable::GetParticleTable()->FindIon(aZ,anN+aZ,0,aZ); |
---|
| 160 | if(aZ == 0 && anN == 1) theDefinition = G4Neutron::Neutron(); |
---|
| 161 | } |
---|
| 162 | else |
---|
| 163 | { |
---|
| 164 | theDefinition = G4ParticleTable::GetParticleTable() |
---|
| 165 | ->FindParticle(output->operator[](particle)->GetPDGCode()); |
---|
| 166 | } |
---|
| 167 | theSec->SetDefinition(theDefinition); |
---|
| 168 | theSec->SetMomentum(output->operator[](particle)->Get4Momentum().vect()); |
---|
| 169 | aResult->AddSecondary(theSec); |
---|
| 170 | delete output->operator[](particle); |
---|
| 171 | } |
---|
| 172 | delete output; |
---|
| 173 | return aResult; |
---|
| 174 | } |
---|
| 175 | |
---|