Changeset 1340 for trunk/source/processes/hadronic/models/util
- Timestamp:
- Nov 5, 2010, 3:45:55 PM (14 years ago)
- Location:
- trunk/source/processes/hadronic/models/util
- Files:
-
- 7 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/source/processes/hadronic/models/util/History
r1315 r1340 1 $Id: History,v 1.40 2010/09/28 17:03:26 vnivanch Exp $ 1 2 ------------------------------------------------------------------- 2 3 … … 15 16 --------------------------------------------------------------- 16 17 18 27 Sep 2010 Vladimir Ivanchenko hadr-mod-util-V09-03-04 19 - G4Fragment - added members numberOfChargedHoles, numberOfShellElectrons 20 and corresponding Get/Set methods; 21 reodered inline methods and extended comments; 22 removed unused private methods and headers 23 24 25 Sep 2010 Michael Kelsey 25 - G4Fragment - Change "setprecision" to "setw" in operator<<, add null 26 pointer check there as well. 27 - History - Add CVS "Id" string at top of file. 28 29 8 Sep 2010 Gunter Folger hadr-mod-util-V09-03-03 30 - G4DecayStrongResonances: cleanup unused #includes 31 - G4Fancy3DNucleus: add integer (A,Z) Init(A,Z) 32 33 19 May 2010 Gabriele Cosmo hadr-mod-util-V09-03-02 34 - G4Fancy3DNucleus: added missing std:: to call to sort() algorithm in code. 35 17 36 19 May 2010 Vladimir Ivanchenko hadr-mod-util-V09-03-01 18 - G4Fragment -minor speedup by adding member and access method19 to GroundStateMass;37 - G4Fragment: minor speedup by adding member and access method 38 to GroundStateMass. 20 39 21 40 10 May 2010 Vladimir Ivanchenko hadr-mod-util-V09-03-00 -
trunk/source/processes/hadronic/models/util/include/G4DecayStrongResonances.hh
r819 r1340 28 28 #define G4DecayStrongResonances_h 1 29 29 30 #include "G4Fancy3DNucleus.hh"31 #include "G4Nucleon.hh"32 #include "G4Nucleus.hh"33 30 #include "G4KineticTrackVector.hh" 34 #include "G4FragmentVector.hh" 35 #include "G4HadFinalState.hh" 36 #include "G4DynamicParticleVector.hh" 37 #include "G4HadTmpUtil.hh" 31 #include "G4ReactionProductVector.hh" 38 32 39 #include <algorithm> 33 class G4V3DNucleus; 40 34 41 35 class G4DecayStrongResonances … … 52 46 53 47 public: 54 G4ReactionProductVector* Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* ) 55 { 56 // decay the strong resonances 57 //static int call_count = 0; 58 //if(call_count++<10) 59 //{ 60 // G4cout << "Security print-out: Entering G4DecayStrongResonances::Propagate"; 61 //} 62 G4ReactionProductVector * theResult; 63 try 64 { 65 theResult = new G4ReactionProductVector; 66 } 67 catch(...) 68 { 69 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: out of memory "); 70 } 71 G4KineticTrackVector *result1, *secondaries, *result; 72 if(!theSecondaries) 73 { 74 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: 0x0 input vector "); 75 } 76 result1=theSecondaries; 77 try 78 { 79 result=new G4KineticTrackVector(); 80 } 81 catch(...) 82 { 83 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: out of memory in "); 84 } 85 86 size_t aResult=0; 87 for (aResult=0; aResult < result1->size(); aResult++) 88 { 89 G4ParticleDefinition * pdef; 90 if(!result1->operator[](aResult)) 91 { 92 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: null pointer in input vector!!! "); 93 } 94 pdef=result1->operator[](aResult)->GetDefinition(); 95 if(!pdef) 96 { 97 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: 0x0 particle definition "); 98 } 99 secondaries=NULL; 100 if ( pdef->GetPDGWidth() > 0 && pdef->GetPDGLifeTime() < 5E-17*s ) 101 { 102 try 103 { 104 secondaries = result1->operator[](aResult)->Decay(); 105 } 106 catch(...) 107 { 108 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: failing in Decay "); 109 } 110 } 111 if ( secondaries == NULL ) 112 { 113 try 114 { 115 result->push_back(result1->operator[](aResult)); 116 } 117 catch(...) 118 { 119 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: push_back failed - out of memory "); 120 } 121 result1->operator[](aResult)=NULL; //protect for clearAndDestroy 122 } 123 else 124 { 125 for (size_t aSecondary=0; aSecondary<secondaries->size(); aSecondary++) 126 { 127 try 128 { 129 result1->push_back(secondaries->operator[](aSecondary)); 130 } 131 catch(...) 132 { 133 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: push_back 1 failed - out of mem"); 134 } 135 } 136 if(secondaries) delete secondaries; 137 } 138 } 139 try 140 { 141 std::for_each(result1->begin(), result1->end(), G4Delete()); 142 delete result1; 143 } 144 catch(...) 145 { 146 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: memory corruption."); 147 } 148 149 // translate to ReactionProducts 150 G4ReactionProduct * it = NULL; 151 for(aResult=0; aResult < result->size(); aResult++) 152 { 153 try 154 { 155 it = new G4ReactionProduct(); 156 } 157 catch(...) 158 { 159 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: out of memory "); 160 } 161 it->SetDefinition((*result)[aResult]->GetDefinition()); 162 it->SetMass((*result)[aResult]->GetDefinition()->GetPDGMass()); 163 it->SetTotalEnergy((*result)[aResult]->Get4Momentum().t()); 164 it->SetMomentum((*result)[aResult]->Get4Momentum().vect()); 165 166 try 167 { 168 theResult->push_back(it); 169 } 170 catch(...) 171 { 172 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: push to result failed - out of mem."); 173 } 174 } 175 try 176 { 177 std::for_each(result->begin(), result->end(), G4Delete()); 178 delete result; 179 } 180 catch(...) 181 { 182 throw G4HadronicException(__FILE__, __LINE__, "DecayStrongRes: memory corruption at end."); 183 } 184 return theResult; 185 } 48 G4ReactionProductVector* Propagate(G4KineticTrackVector* theSecondaries, 49 G4V3DNucleus* ); 186 50 }; 187 51 -
trunk/source/processes/hadronic/models/util/include/G4Fancy3DNucleus.hh
r1196 r1340 45 45 #include <vector> 46 46 47 // to test if we can drop old interface for (A,Z), comment next line.. 48 //#define NON_INTEGER_A_Z 1 49 47 50 class G4Fancy3DNucleus : public G4V3DNucleus 48 51 { … … 67 70 68 71 public: 72 #if defined(NON_INTEGER_A_Z) 69 73 void Init(G4double theA, G4double theZ); 74 #endif 75 void Init(G4int theA, G4int theZ); 70 76 G4bool StartLoop(); 71 77 G4Nucleon * GetNextNucleon(); -
trunk/source/processes/hadronic/models/util/include/G4FermiMomentum.hh
r819 r1340 40 40 ~G4FermiMomentum(); 41 41 42 inline void Init(G4 double anA, G4doubleaZ) {theA = anA; theZ = aZ;}42 inline void Init(G4int anA, G4int aZ) {theA = anA; theZ = aZ;} 43 43 44 44 inline G4double GetFermiMomentum(G4double density) … … 67 67 private: 68 68 69 G4 doubletheA;70 G4 doubletheZ;69 G4int theA; 70 G4int theZ; 71 71 // pmax= hbar * c * ( 3* pi**2 * rho )**(1/3) = 72 72 // hbar * c * ( 3* pi**2 )**(1/3) * rho**(1/3)= -
trunk/source/processes/hadronic/models/util/include/G4Fragment.hh
r1337 r1340 24 24 // ******************************************************************** 25 25 // 26 // $Id: G4Fragment.hh,v 1.1 1 2010/05/19 10:23:00 vnivanch Exp $27 // GEANT4 tag $Name: geant4-09-0 4-beta-01$26 // $Id: G4Fragment.hh,v 1.16 2010/09/28 16:09:00 vnivanch Exp $ 27 // GEANT4 tag $Name: geant4-09-03-ref-09 $ 28 28 // 29 29 //--------------------------------------------------------------------- … … 41 41 // method which allowing to compute this value once and use 42 42 // many times 43 // 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes 44 // and neutron holes as members of the class and Get/Set methods; 45 // removed not needed 'const'; removed old debug staff and unused 46 // private methods; add comments and reorder methods for 47 // better reading 43 48 44 49 #ifndef G4Fragment_h 45 50 #define G4Fragment_h 1 46 51 47 #include "G4ios.hh"48 #include <iomanip>49 52 #include <vector> 50 53 51 54 #include "globals.hh" 52 55 #include "G4LorentzVector.hh" 53 //#include "G4ParticleMomentum.hh"54 56 #include "G4ThreeVector.hh" 55 57 #include "G4NucleiProperties.hh" 56 #include "G4ParticleTable.hh"57 #include "G4IonTable.hh"58 //#include "G4ParticleTable.hh" 59 //#include "G4IonTable.hh" 58 60 #include "Randomize.hh" 59 61 #include "G4Proton.hh" 60 62 #include "G4Neutron.hh" 61 #include "G4HadronicException.hh"62 63 #include "G4HadTmpUtil.hh" 63 64 64 65 65 class G4ParticleDefinition; 66 66 67 class G4Fragment; // Forward deckaration67 class G4Fragment; 68 68 typedef std::vector<G4Fragment*> G4FragmentVector; 69 69 … … 84 84 85 85 // A,Z and 4-momentum - main constructor for fragment 86 G4Fragment(const G4int A, const G4int Z, const G4LorentzVector& aMomentum); 87 88 // 4-momentum and pointer to G4particleDefinition (for gammas) 89 G4Fragment(const G4LorentzVector& aMomentum, G4ParticleDefinition * aParticleDefinition); 86 G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum); 87 88 // 4-momentum and pointer to G4particleDefinition (for gammas, e-) 89 G4Fragment(const G4LorentzVector& aMomentum, 90 G4ParticleDefinition* aParticleDefinition); 90 91 91 92 // ============= OPERATORS ================== … … 98 99 friend std::ostream& operator<<(std::ostream&, const G4Fragment&); 99 100 100 // ============= METHODS ================== 101 102 inline G4double GetA() const; 103 inline void SetA(const G4double value); 104 105 inline G4double GetZ() const; 106 inline void SetZ(const G4double value); 101 // ============= GENERAL METHODS ================== 107 102 108 103 inline G4int GetZ_asInt() const; … … 111 106 112 107 inline G4double GetExcitationEnergy() const; 113 void SetExcitationEnergy(const G4double value);114 115 inline const G4LorentzVector& GetMomentum() const;116 inline void SetMomentum(const G4LorentzVector& value);117 118 inline const G4ThreeVector& GetAngularMomentum() const;119 inline void SetAngularMomentum(const G4ThreeVector& value);120 121 inline G4int GetNumberOfExcitons() const;122 123 inline G4int GetNumberOfHoles() const;124 inline void SetNumberOfHoles(const G4int value);125 126 inline G4int GetNumberOfCharged() const;127 void SetNumberOfCharged(const G4int value);128 129 inline G4int GetNumberOfParticles() const;130 inline void SetNumberOfParticles(const G4int value);131 132 inline G4ParticleDefinition * GetParticleDefinition() const;133 inline void SetParticleDefinition(G4ParticleDefinition * aParticleDefinition);134 135 inline G4double GetCreationTime() const;136 inline void SetCreationTime(const G4double time);137 108 138 109 inline G4double GetGroundStateMass() const; 139 110 140 111 inline G4double GetBindingEnergy() const; 112 113 inline const G4LorentzVector& GetMomentum() const; 114 inline void SetMomentum(const G4LorentzVector& value); 115 116 inline const G4ThreeVector& GetAngularMomentum() const; 117 inline void SetAngularMomentum(const G4ThreeVector& value); 141 118 142 119 // computation of mass for any Z and A 143 inline G4double ComputeGroundStateMass(const G4int Z, const G4int A) const; 144 145 #ifdef PRECOMPOUND_TEST 146 G4String GetCreatorModel() const; 147 void SetCreatorModel(const G4String & aModel); 148 #endif 120 inline G4double ComputeGroundStateMass(G4int Z, G4int A) const; 121 122 // obsolete methods 123 inline G4double GetZ() const; 124 inline G4double GetA() const; 125 inline void SetZ(G4double value); 126 inline void SetA(G4double value); 127 128 // ============= METHODS FOR PRE-COMPOUND MODEL =============== 129 130 inline G4int GetNumberOfExcitons() const; 131 132 inline G4int GetNumberOfParticles() const; 133 inline G4int GetNumberOfCharged() const; 134 inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP); 135 136 inline G4int GetNumberOfHoles() const; 137 inline G4int GetNumberOfChargedHoles() const; 138 inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0); 139 140 // these methods will be removed in future 141 inline void SetNumberOfParticles(G4int value); 142 inline void SetNumberOfCharged(G4int value); 143 144 // ============= METHODS FOR PHOTON EVAPORATION =============== 145 146 inline G4int GetNumberOfElectrons() const; 147 inline void SetNumberOfElectrons(G4int value); 148 149 inline G4ParticleDefinition * GetParticleDefinition() const; 150 inline void SetParticleDefinition(G4ParticleDefinition * p); 151 152 inline G4double GetCreationTime() const; 153 inline void SetCreationTime(G4double time); 154 155 // ============= PRIVATE METHODS ============================== 149 156 150 157 private: 151 158 152 void ExcitationEnegryWarning(); 159 void ExcitationEnergyWarning(); 160 161 void NumberOfExitationWarning(const G4String&); 162 163 inline void CalculateExcitationEnergy(); 153 164 154 165 inline void CalculateGroundStateMass(); 155 166 156 inline void CalculateExcitationEnergy();157 158 G4ThreeVector IsotropicRandom3Vector(const G4double Magnitude = 1.0) const;159 160 167 // ============= DATA MEMBERS ================== 161 168 … … 173 180 174 181 G4ThreeVector theAngularMomentum; 182 183 // Exciton model data members 175 184 176 185 G4int numberOfParticles; 177 186 187 G4int numberOfCharged; 188 178 189 G4int numberOfHoles; 179 190 180 G4int numberOfCharged; 181 182 // Gamma evaporation requeriments 191 G4int numberOfChargedHoles; 192 193 // Gamma evaporation data members 194 195 G4int numberOfShellElectrons; 183 196 184 197 G4ParticleDefinition * theParticleDefinition; … … 186 199 G4double theCreationTime; 187 200 188 #ifdef PRECOMPOUND_TEST189 G4String theCreatorModel;190 #endif191 201 }; 192 202 193 // Class G4Fragment 203 // ============= INLINE METHOD IMPLEMENTATIONS =================== 204 205 inline void G4Fragment::CalculateExcitationEnergy() 206 { 207 theExcitationEnergy = theMomentum.mag() - theGroundStateMass; 208 if(theExcitationEnergy < 0.0) { ExcitationEnergyWarning(); } 209 } 210 194 211 inline void G4Fragment::CalculateGroundStateMass() 195 212 { 196 213 theGroundStateMass = G4NucleiProperties::GetNuclearMass(theA, theZ); 197 }198 199 inline G4double G4Fragment::GetA() const200 {201 return G4double(theA);202 }203 204 inline void G4Fragment::SetA(const G4double value)205 {206 theA = G4lrint(value);207 CalculateGroundStateMass();208 }209 210 inline G4double G4Fragment::GetZ() const211 {212 return G4double(theZ);213 }214 215 inline void G4Fragment::SetZ(const G4double value)216 {217 theZ = G4lrint(value);218 CalculateGroundStateMass();219 214 } 220 215 … … 238 233 inline G4double G4Fragment::GetExcitationEnergy() const 239 234 { 240 // temporary fix for what seems to be241 // a problem with rounding errors for on-shell lorentz-vectors in CLHEP.242 // HPW Apr 1999 @@@@@@@243 244 //VI if(std::abs(theExcitationEnergy)<10*eV) return 0;245 235 return theExcitationEnergy; 246 236 } 247 237 248 inline const G4LorentzVector& G4Fragment::GetMomentum() const249 {250 return theMomentum;251 }252 253 inline const G4ThreeVector& G4Fragment::GetAngularMomentum() const254 {255 return theAngularMomentum;256 }257 258 inline void G4Fragment::SetAngularMomentum(const G4ThreeVector& value)259 {260 theAngularMomentum = value;261 }262 263 inline G4int G4Fragment::GetNumberOfExcitons() const264 {265 return numberOfParticles + numberOfHoles;266 }267 268 inline void G4Fragment::SetNumberOfParticles(const G4int value)269 {270 numberOfParticles = value;271 }272 273 inline G4int G4Fragment::GetNumberOfHoles() const274 {275 return numberOfHoles;276 }277 278 inline void G4Fragment::SetNumberOfHoles(const G4int value)279 {280 numberOfHoles = value;281 }282 283 inline G4int G4Fragment::GetNumberOfCharged() const284 {285 return numberOfCharged;286 }287 288 inline void G4Fragment::SetNumberOfCharged(const G4int value)289 {290 if (value <= numberOfParticles) { numberOfCharged = value; }291 else292 {293 G4String text = "G4Fragment::SetNumberOfCharged: Number of charged particles can't be greater than number of particles";294 throw G4HadronicException(__FILE__, __LINE__, text);295 }296 }297 298 inline G4int G4Fragment::GetNumberOfParticles() const299 {300 return numberOfParticles;301 }302 303 inline G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const304 {305 return theParticleDefinition;306 }307 308 inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * aParticleDefinition)309 {310 theParticleDefinition = aParticleDefinition;311 }312 313 inline G4double G4Fragment::GetCreationTime() const314 {315 return theCreationTime;316 }317 318 inline void G4Fragment::SetCreationTime(const G4double time)319 {320 theCreationTime = time;321 }322 323 238 inline G4double G4Fragment::GetGroundStateMass() const 324 239 { … … 326 241 } 327 242 328 inline G4double329 G4Fragment::ComputeGroundStateMass(const G4int Z, const G4int A) const330 {331 return G4NucleiProperties::GetNuclearMass(A, Z);332 }333 334 inline void G4Fragment::CalculateExcitationEnergy()335 {336 theExcitationEnergy = theMomentum.mag() - theGroundStateMass;337 if(theExcitationEnergy < 0.0) { ExcitationEnegryWarning(); }338 }339 340 243 inline G4double G4Fragment::GetBindingEnergy() const 341 244 { … … 344 247 } 345 248 249 inline const G4LorentzVector& G4Fragment::GetMomentum() const 250 { 251 return theMomentum; 252 } 253 346 254 inline void G4Fragment::SetMomentum(const G4LorentzVector& value) 347 255 { … … 350 258 } 351 259 260 inline const G4ThreeVector& G4Fragment::GetAngularMomentum() const 261 { 262 return theAngularMomentum; 263 } 264 265 inline void G4Fragment::SetAngularMomentum(const G4ThreeVector& value) 266 { 267 theAngularMomentum = value; 268 } 269 270 inline G4double 271 G4Fragment::ComputeGroundStateMass(G4int Z, G4int A) const 272 { 273 return G4NucleiProperties::GetNuclearMass(A, Z); 274 } 275 276 inline G4double G4Fragment::GetZ() const 277 { 278 return G4double(theZ); 279 } 280 281 inline G4double G4Fragment::GetA() const 282 { 283 return G4double(theA); 284 } 285 286 inline void G4Fragment::SetZ(const G4double value) 287 { 288 theZ = G4lrint(value); 289 CalculateGroundStateMass(); 290 } 291 292 inline void G4Fragment::SetA(const G4double value) 293 { 294 theA = G4lrint(value); 295 CalculateGroundStateMass(); 296 } 297 298 inline G4int G4Fragment::GetNumberOfExcitons() const 299 { 300 return numberOfParticles + numberOfHoles; 301 } 302 303 inline G4int G4Fragment::GetNumberOfParticles() const 304 { 305 return numberOfParticles; 306 } 307 308 inline G4int G4Fragment::GetNumberOfCharged() const 309 { 310 return numberOfCharged; 311 } 312 313 inline 314 void G4Fragment::SetNumberOfExcitedParticle(G4int valueTot, G4int valueP) 315 { 316 numberOfParticles = valueTot; 317 numberOfCharged = valueP; 318 if(valueTot < valueP) { 319 NumberOfExitationWarning("SetNumberOfExcitedParticle"); 320 } 321 } 322 323 inline G4int G4Fragment::GetNumberOfHoles() const 324 { 325 return numberOfHoles; 326 } 327 328 inline G4int G4Fragment::GetNumberOfChargedHoles() const 329 { 330 return numberOfChargedHoles; 331 } 332 333 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP) 334 { 335 numberOfHoles = valueTot; 336 numberOfChargedHoles = valueP; 337 if(valueTot < valueP) { 338 NumberOfExitationWarning("SetNumberOfHoles"); 339 } 340 } 341 342 inline void G4Fragment::SetNumberOfParticles(G4int value) 343 { 344 numberOfParticles = value; 345 } 346 347 inline void G4Fragment::SetNumberOfCharged(G4int value) 348 { 349 numberOfCharged = value; 350 if(value > numberOfParticles) { 351 NumberOfExitationWarning("SetNumberOfCharged"); 352 } 353 } 354 355 inline G4int G4Fragment::GetNumberOfElectrons() const 356 { 357 return numberOfShellElectrons; 358 } 359 360 inline void G4Fragment::SetNumberOfElectrons(G4int value) 361 { 362 numberOfShellElectrons = value; 363 } 364 365 inline 366 G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const 367 { 368 return theParticleDefinition; 369 } 370 371 inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * p) 372 { 373 theParticleDefinition = p; 374 } 375 376 inline G4double G4Fragment::GetCreationTime() const 377 { 378 return theCreationTime; 379 } 380 381 inline void G4Fragment::SetCreationTime(G4double time) 382 { 383 theCreationTime = time; 384 } 385 352 386 #endif 353 387 -
trunk/source/processes/hadronic/models/util/src/G4Fancy3DNucleus.cc
r1196 r1340 61 61 } 62 62 63 63 #if defined(NON_INTEGER_A_Z) 64 64 void G4Fancy3DNucleus::Init(G4double theA, G4double theZ) 65 { 66 G4int intZ = G4int(theZ); 67 G4int intA= ( G4UniformRand()>theA-G4int(theA) ) ? G4int(theA) : G4int(theA)+1; 68 // forward to integer Init() 69 Init(intA, intZ); 70 71 } 72 #endif 73 74 void G4Fancy3DNucleus::Init(G4int theA, G4int theZ) 65 75 { 66 76 // G4cout << "G4Fancy3DNucleus::Init(theA, theZ) called"<<G4endl; … … 70 80 theRWNucleons.clear(); 71 81 72 myZ = G4int(theZ);73 myA= ( G4UniformRand()>theA-G4int(theA) ) ? G4int(theA) : G4int(theA)+1;82 myZ = theZ; 83 myA= theA; 74 84 75 85 theNucleons = new G4Nucleon[myA]; … … 144 154 if (theRWNucleons.size() < 2 ) return; 145 155 146 s ort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ);156 std::sort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ); 147 157 148 158 // now copy sorted nucleons to theNucleons array. TheRWNucleons are pointers in theNucleons … … 169 179 170 180 if (theRWNucleons.size() < 2 ) return; 171 s ort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ);181 std::sort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ); 172 182 173 183 // now copy sorted nucleons to theNucleons array. TheRWNucleons are pointers in theNucleons -
trunk/source/processes/hadronic/models/util/src/G4Fragment.cc
r1337 r1340 24 24 // ******************************************************************** 25 25 // 26 // $Id: G4Fragment.cc,v 1. 16 2010/05/18 18:52:07vnivanch Exp $27 // GEANT4 tag $Name: geant4-09-0 4-beta-01$26 // $Id: G4Fragment.cc,v 1.21 2010/09/28 16:06:32 vnivanch Exp $ 27 // GEANT4 tag $Name: geant4-09-03-ref-09 $ 28 28 // 29 29 //--------------------------------------------------------------------- … … 37 37 // 03.05.2010 V.Ivanchenko General cleanup; moved obsolete methods from 38 38 // inline to source 39 // 39 // 25.09.2010 M. Kelsey -- Change "setprecision" to "setwidth" in printout, 40 // add null pointer check. 40 41 41 42 #include "G4Fragment.hh" 42 43 #include "G4HadronicException.hh" 43 #include "G4HadTmpUtil.hh"44 44 #include "G4Gamma.hh" 45 45 #include "G4Electron.hh" 46 #include "G4ios.hh" 47 #include <iomanip> 46 48 47 49 G4int G4Fragment::errCount = 0; … … 53 55 theExcitationEnergy(0.0), 54 56 theGroundStateMass(0.0), 55 theMomentum( 0),56 theAngularMomentum( 0),57 theMomentum(G4LorentzVector(0,0,0,0)), 58 theAngularMomentum(G4ThreeVector(0,0,0)), 57 59 numberOfParticles(0), 60 numberOfCharged(0), 58 61 numberOfHoles(0), 59 numberOfCharged(0), 62 numberOfChargedHoles(0), 63 numberOfShellElectrons(0), 60 64 theParticleDefinition(0), 61 65 theCreationTime(0.0) 62 #ifdef PRECOMPOUND_TEST63 ,theCreatorModel("No name")64 #endif65 66 {} 66 67 … … 75 76 theAngularMomentum = right.theAngularMomentum; 76 77 numberOfParticles = right.numberOfParticles; 78 numberOfCharged = right.numberOfCharged; 77 79 numberOfHoles = right.numberOfHoles; 78 numberOfCharged = right.numberOfCharged; 80 numberOfChargedHoles = right.numberOfChargedHoles; 81 numberOfShellElectrons = right.numberOfShellElectrons; 79 82 theParticleDefinition = right.theParticleDefinition; 80 83 theCreationTime = right.theCreationTime; 81 #ifdef PRECOMPOUND_TEST82 theCreatorModel = right.theCreatorModel;83 #endif84 84 } 85 85 … … 87 87 {} 88 88 89 G4Fragment::G4Fragment( const G4int A, constG4int Z, const G4LorentzVector& aMomentum) :89 G4Fragment::G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum) : 90 90 theA(A), 91 91 theZ(Z), 92 92 theMomentum(aMomentum), 93 theAngularMomentum( 0),93 theAngularMomentum(G4ThreeVector(0,0,0)), 94 94 numberOfParticles(0), 95 numberOfCharged(0), 95 96 numberOfHoles(0), 96 numberOfCharged(0), 97 numberOfChargedHoles(0), 98 numberOfShellElectrons(0), 97 99 theParticleDefinition(0), 98 100 theCreationTime(0.0) 99 #ifdef PRECOMPOUND_TEST100 ,theCreatorModel("No name")101 #endif102 101 { 103 102 theExcitationEnergy = 0.0; … … 107 106 CalculateExcitationEnergy(); 108 107 } 109 /* 110 theExcitationEnergy = theMomentum.mag() - 111 G4ParticleTable::GetParticleTable()->GetIonTable() 112 ->GetIonMass( G4lrint(theZ), G4lrint(theA) ); 113 if (theExcitationEnergy < 0.0) { 114 if (theExcitationEnergy > -10.0 * eV || 0 == G4lrint(theA)) { 115 theExcitationEnergy = 0.0; 116 } else { 117 G4cout << "A, Z, momentum, theExcitationEnergy"<< 118 A<<" "<<Z<<" "<<aMomentum<<" "<<theExcitationEnergy<<G4endl; 119 G4String text = "G4Fragment::G4Fragment Excitation Energy < 0.0!"; 120 throw G4HadronicException(__FILE__, __LINE__, text); 121 } 122 } 123 */ 124 } 125 108 } 126 109 127 110 // This constructor is for initialize photons or electrons … … 131 114 theZ(0), 132 115 theMomentum(aMomentum), 133 theAngularMomentum( 0),116 theAngularMomentum(G4ThreeVector(0,0,0)), 134 117 numberOfParticles(0), 118 numberOfCharged(0), 135 119 numberOfHoles(0), 136 numberOfCharged(0), 120 numberOfChargedHoles(0), 121 numberOfShellElectrons(0), 137 122 theParticleDefinition(aParticleDefinition), 138 123 theCreationTime(0.0) 139 #ifdef PRECOMPOUND_TEST140 ,theCreatorModel("No name")141 #endif142 124 { 143 125 theExcitationEnergy = 0.0; … … 161 143 theAngularMomentum = right.theAngularMomentum; 162 144 numberOfParticles = right.numberOfParticles; 145 numberOfCharged = right.numberOfCharged; 163 146 numberOfHoles = right.numberOfHoles; 164 numberOfCharged = right.numberOfCharged; 147 numberOfChargedHoles = right.numberOfChargedHoles; 148 numberOfShellElectrons = right.numberOfShellElectrons; 165 149 theParticleDefinition = right.theParticleDefinition; 166 150 theCreationTime = right.theCreationTime; 167 #ifdef PRECOMPOUND_TEST168 theCreatorModel = right.theCreatorModel;169 #endif170 151 } 171 152 return *this; … … 184 165 std::ostream& operator << (std::ostream &out, const G4Fragment *theFragment) 185 166 { 167 if (!theFragment) { 168 out << "Fragment: null pointer "; 169 return out; 170 } 171 186 172 std::ios::fmtflags old_floatfield = out.flags(); 187 173 out.setf(std::ios::floatfield); 188 174 189 out 190 << "Fragment: A = " << std::setprecision(3) << theFragment->theA 191 << ", Z = " << std::setprecision(3) << theFragment->theZ ; 175 out << "Fragment: A = " << std::setw(3) << theFragment->theA 176 << ", Z = " << std::setw(3) << theFragment->theZ ; 192 177 out.setf(std::ios::scientific,std::ios::floatfield); 193 out 194 << ", U = " << theFragment->GetExcitationEnergy()/MeV 195 << " MeV" << G4endl 196 << " P = (" 197 << theFragment->theMomentum.x()/MeV << "," 198 << theFragment->theMomentum.y()/MeV << "," 199 << theFragment->theMomentum.z()/MeV 200 << ") MeV E = " 201 << theFragment->theMomentum.t()/MeV << " MeV"; 202 178 179 // Store user's precision setting and reset to (3) here: back-compatibility 180 std::streamsize floatPrec = out.precision(); 181 182 out << std::setprecision(3) 183 << ", U = " << theFragment->GetExcitationEnergy()/CLHEP::MeV 184 << " MeV" << G4endl 185 << " P = (" 186 << theFragment->theMomentum.x()/CLHEP::MeV << "," 187 << theFragment->theMomentum.y()/CLHEP::MeV << "," 188 << theFragment->theMomentum.z()/CLHEP::MeV 189 << ") MeV E = " 190 << theFragment->theMomentum.t()/CLHEP::MeV << " MeV" 191 << G4endl; 192 203 193 // What about Angular momentum??? 204 194 205 195 if (theFragment->GetNumberOfExcitons() != 0) { 206 out << G4endl;207 196 out << " " 208 << "#Particles = " << theFragment->numberOfParticles 209 << ", #Holes = " << theFragment->numberOfHoles 210 << ", #Charged = " << theFragment->numberOfCharged; 197 << "#Particles= " << theFragment->numberOfParticles 198 << ", #Charged= " << theFragment->numberOfCharged 199 << ", #Holes= " << theFragment->numberOfHoles 200 << ", #ChargedHoles= " << theFragment->numberOfChargedHoles 201 << G4endl; 211 202 } 212 203 out.setf(old_floatfield,std::ios::floatfield); 204 out.precision(floatPrec); 213 205 214 206 return out; 215 216 207 } 217 208 … … 222 213 } 223 214 224 void G4Fragment::ExcitationEne gryWarning()225 { 226 if (theExcitationEnergy < -10 .0 *eV) {215 void G4Fragment::ExcitationEnergyWarning() 216 { 217 if (theExcitationEnergy < -10 * CLHEP::eV) { 227 218 ++errCount; 228 219 if ( errCount <= 10 ) { 229 220 G4cout << "G4Fragment::CalculateExcitationEnergy(): Excitation Energy = " 230 << theExcitationEnergy/ MeV << " MeV for A = "231 << theA << " and Z= " << theZ << G4endl;221 << theExcitationEnergy/CLHEP::MeV << " MeV for A = " 222 << theA << " and Z= " << theZ << G4endl; 232 223 if( errCount == 10 ) { 233 224 G4String text = "G4Fragment::G4Fragment Excitation Energy < 0.0!"; … … 239 230 } 240 231 241 G4ThreeVector G4Fragment::IsotropicRandom3Vector(const G4double Magnitude) const 242 // Create a unit vector with a random direction isotropically distributed 243 { 244 G4double CosTheta = 1.0 - 2.0*G4UniformRand(); 245 G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta); 246 G4double Phi = twopi*G4UniformRand(); 247 G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta, 248 Magnitude*std::sin(Phi)*SinTheta, 249 Magnitude*CosTheta); 250 251 return Vector; 252 } 253 254 void G4Fragment::SetExcitationEnergy(const G4double ) 255 { 256 // theExcitationEnergy = value; 257 G4cout << "Warning: G4Fragment::SetExcitationEnergy() is a dummy method. Please, avoid to use it." << G4endl; 258 } 259 260 #ifdef PRECOMPOUND_TEST 261 G4String G4Fragment::GetCreatorModel() const 262 { 263 return theCreatorModel; 264 } 265 266 void G4Fragment::SetCreatorModel(const G4String & aModel) 267 { 268 theCreatorModel = aModel; 269 } 270 #endif 232 void G4Fragment::NumberOfExitationWarning(const G4String& value) 233 { 234 G4cout << "G4Fragment::"<< value << " ERROR " 235 << G4endl; 236 G4cout << this << G4endl; 237 G4String text = "G4Fragment::G4Fragment wrong exciton number "; 238 throw G4HadronicException(__FILE__, __LINE__, text); 239 }
Note: See TracChangeset
for help on using the changeset viewer.