Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:18

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 //---------------------------------------------------------------------
0028 //
0029 // Geant4 header G4Fragment
0030 //
0031 // by V. Lara (May 1998)
0032 //
0033 // Modifications:
0034 // 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects 
0035 //            are accessed by reference; remove double return 
0036 //            tolerance of excitation energy at modent it is computed;
0037 //            safe computation of excitation for exotic fragments
0038 // 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline
0039 //            method which allowing to compute this value once and use 
0040 //            many times
0041 // 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes
0042 //            and neutron holes as members of the class and Get/Set methods;
0043 //            removed not needed 'const'; removed old debug staff and unused
0044 //            private methods; add comments and reorder methods for 
0045 //            better reading
0046 // 27.10.2021 A.Ribon extension for hypernuclei.
0047 
0048 #ifndef G4Fragment_h
0049 #define G4Fragment_h 1
0050 
0051 #include "globals.hh"
0052 #include "G4Allocator.hh"
0053 #include "G4LorentzVector.hh"
0054 #include "G4ThreeVector.hh"
0055 #include "G4NuclearPolarization.hh"
0056 #include "G4NucleiProperties.hh"
0057 #include "G4HyperNucleiProperties.hh"
0058 #include "G4Proton.hh"
0059 #include "G4Neutron.hh"
0060 #include <vector>
0061 
0062 class G4ParticleDefinition;
0063 
0064 class G4Fragment;     
0065 typedef std::vector<G4Fragment*> G4FragmentVector;
0066 
0067 class G4Fragment 
0068 {
0069 public:
0070 
0071   // ============= CONSTRUCTORS ==================
0072 
0073   // Default constructor - obsolete
0074   G4Fragment();
0075 
0076   // Destructor
0077   ~G4Fragment();
0078 
0079   // Copy constructor
0080   G4Fragment(const G4Fragment &right);
0081 
0082   // A,Z and 4-momentum - main constructor for fragment
0083   G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum);
0084 
0085   // A,Z,numberOfLambdas and 4-momentum
0086   G4Fragment(G4int A, G4int Z, G4int numberOfLambdas,
0087              const G4LorentzVector& aMomentum);
0088 
0089   // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
0090   G4Fragment(const G4LorentzVector& aMomentum, 
0091          const G4ParticleDefinition* aParticleDefinition);
0092 
0093   // ============= OPERATORS ==================
0094     
0095   G4Fragment & operator=(const G4Fragment &right);
0096   G4bool operator==(const G4Fragment &right) const;
0097   G4bool operator!=(const G4Fragment &right) const;
0098 
0099   friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
0100 
0101   //  new/delete operators are overloded to use G4Allocator
0102   inline void *operator new(size_t);
0103   inline void operator delete(void *aFragment);
0104 
0105   // ============= GENERAL METHODS ==================
0106 
0107   inline G4int GetZ_asInt() const;
0108   inline G4int GetA_asInt() const;
0109 
0110   // update number of nucleons without check on input
0111   // ground state mass is not recomputed
0112   inline void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0);
0113 
0114   // non-negative number of lambdas/anti-lambdas
0115   // ground state mass is not recomputed
0116   void SetNumberOfLambdas(G4int numberOfLambdas);
0117 
0118   inline G4int GetNumberOfLambdas() const;
0119 
0120   inline G4double GetExcitationEnergy() const;
0121 
0122   inline G4double GetGroundStateMass() const;
0123    
0124   inline const G4LorentzVector& GetMomentum() const;
0125 
0126   // ground state mass and excitation energy are recomputed
0127   inline G4double RecomputeGroundStateMass();
0128 
0129   // update main fragment parameters full check on input 
0130   // ground state mass and excitation energy are recomputed
0131   inline void SetMomentum(const G4LorentzVector& value);
0132   inline void SetZAandMomentum(const G4LorentzVector&,
0133                                G4int Z, G4int A,
0134                                G4int nLambdas = 0);
0135 
0136   // ground state mass is not recomputed
0137   void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector&);
0138   G4double GetBindingEnergy() const;
0139     
0140   // computation of mass for any imput Z, A and number of Lambdas
0141   // no check on input values
0142   inline G4double ComputeGroundStateMass(G4int Z, G4int A,
0143                                          G4int nLambdas = 0) const;
0144 
0145   // extra methods
0146   inline G4double GetSpin() const;
0147   inline void SetSpin(G4double value);
0148 
0149   inline G4int GetCreatorModelID() const;
0150   inline void SetCreatorModelID(G4int value);
0151 
0152   inline G4bool IsLongLived() const;
0153   inline void SetLongLived(G4bool value);
0154 
0155   // obsolete methods
0156   inline G4double GetZ() const;
0157   inline G4double GetA() const;
0158   inline void SetZ(G4double value);
0159   inline void SetA(G4double value);
0160   
0161   // ============= METHODS FOR PRE-COMPOUND MODEL ===============
0162 
0163   inline G4int GetNumberOfExcitons() const;
0164   
0165   inline G4int GetNumberOfParticles() const;
0166   inline G4int GetNumberOfCharged() const;
0167   inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
0168 
0169   inline G4int GetNumberOfHoles() const;
0170   inline G4int GetNumberOfChargedHoles() const;
0171   inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
0172   
0173   // these methods will be removed in future
0174   inline void SetNumberOfParticles(G4int value);
0175   inline void SetNumberOfCharged(G4int value);
0176 
0177   // ============= METHODS FOR PHOTON EVAPORATION ===============
0178 
0179   inline G4int GetNumberOfElectrons() const;
0180   inline void SetNumberOfElectrons(G4int value);
0181 
0182   inline G4int GetFloatingLevelNumber() const;
0183   inline void SetFloatingLevelNumber(G4int value);
0184 
0185   inline const G4ParticleDefinition * GetParticleDefinition() const;
0186   inline void SetParticleDefinition(const G4ParticleDefinition * p);
0187 
0188   inline G4double GetCreationTime() const;
0189   inline void SetCreationTime(G4double time);
0190 
0191   // G4Fragment class is not responsible for creation and delition of 
0192   // G4NuclearPolarization object
0193   inline G4NuclearPolarization* NuclearPolarization();
0194   inline G4NuclearPolarization* GetNuclearPolarization() const;
0195   inline void SetNuclearPolarization(G4NuclearPolarization*);
0196 
0197   void SetAngularMomentum(const G4ThreeVector&);
0198   G4ThreeVector GetAngularMomentum() const;
0199 
0200   // ============= PRIVATE METHODS ==============================
0201 
0202 private:
0203 
0204   void CalculateMassAndExcitationEnergy();
0205 
0206   void ExcitationEnergyWarning();
0207 
0208   void NumberOfExitationWarning(const G4String&);
0209 
0210   // ============= DATA MEMBERS ==================
0211 
0212   G4int theA;
0213   
0214   G4int theZ;
0215 
0216   // Non-negative number of lambdas/anti-lambdas inside the nucleus/anti-nucleus
0217   G4int theL;  
0218   
0219   G4double theExcitationEnergy;
0220 
0221   G4double theGroundStateMass;
0222 
0223   G4LorentzVector theMomentum;
0224   
0225   // Nuclear polarisation by default is nullptr
0226   G4NuclearPolarization* thePolarization;
0227 
0228   // creator model type
0229   G4int creatorModel;
0230 
0231   // Exciton model data members  
0232   G4int numberOfParticles;  
0233   G4int numberOfCharged;
0234   G4int numberOfHoles;
0235   G4int numberOfChargedHoles;
0236 
0237   // Gamma evaporation data members
0238   G4int numberOfShellElectrons;
0239   G4int xLevel;
0240 
0241   const G4ParticleDefinition* theParticleDefinition;
0242   
0243   G4double spin;
0244   G4double theCreationTime;
0245 
0246   G4bool isLongLived = false; 
0247 };
0248 
0249 // ============= INLINE METHOD IMPLEMENTATIONS ===================
0250 
0251 #if defined G4HADRONIC_ALLOC_EXPORT
0252   extern G4DLLEXPORT G4Allocator<G4Fragment>*& pFragmentAllocator();
0253 #else
0254   extern G4DLLIMPORT G4Allocator<G4Fragment>*& pFragmentAllocator();
0255 #endif
0256 
0257 inline void * G4Fragment::operator new(size_t)
0258 {
0259   if (!pFragmentAllocator()) { 
0260     pFragmentAllocator() = new G4Allocator<G4Fragment>;
0261   }
0262   return (void*) pFragmentAllocator()->MallocSingle();
0263 }
0264 
0265 inline void G4Fragment::operator delete(void * aFragment)
0266 {
0267   pFragmentAllocator()->FreeSingle((G4Fragment *) aFragment);
0268 }
0269 
0270 inline G4double 
0271 G4Fragment::ComputeGroundStateMass(G4int Z, G4int A, G4int nLambdas) const
0272 {
0273   return ( nLambdas <= 0 ) 
0274     ? G4NucleiProperties::GetNuclearMass(A, Z) 
0275     : G4HyperNucleiProperties::GetNuclearMass(A, Z, nLambdas); 
0276 }
0277 
0278 inline G4double G4Fragment::RecomputeGroundStateMass()
0279 {
0280   CalculateMassAndExcitationEnergy();
0281   return theGroundStateMass;
0282 }
0283 
0284 inline G4int G4Fragment::GetA_asInt() const
0285 {
0286   return theA;
0287 }
0288 
0289 inline G4int G4Fragment::GetZ_asInt()  const
0290 {
0291   return theZ;
0292 }
0293 
0294 inline void 
0295 G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew)
0296 {
0297   theZ = Znew;
0298   theA = Anew;
0299   theL = std::max(Lnew, 0);
0300 }
0301 
0302 inline void G4Fragment::SetNumberOfLambdas(G4int Lnew)
0303 {
0304   theL = std::max(Lnew, 0);
0305 }
0306 
0307 inline G4int G4Fragment::GetNumberOfLambdas() const
0308 {
0309   return theL;
0310 }
0311 
0312 inline G4double G4Fragment::GetExcitationEnergy()  const
0313 {
0314   return theExcitationEnergy;
0315 }
0316 
0317 inline G4double G4Fragment::GetGroundStateMass() const
0318 {
0319   return theGroundStateMass; 
0320 }
0321 
0322 inline const G4LorentzVector& G4Fragment::GetMomentum()  const
0323 {
0324   return theMomentum;
0325 }
0326 
0327 inline void G4Fragment::SetMomentum(const G4LorentzVector& value)
0328 {
0329   theMomentum = value;
0330   CalculateMassAndExcitationEnergy();
0331 }
0332 
0333 inline void 
0334 G4Fragment::SetZAandMomentum(const G4LorentzVector& v,
0335                              G4int Z, G4int A, G4int nLambdas)
0336 {
0337   SetZandA_asInt(Z, A, nLambdas);
0338   SetMomentum(v);
0339 }
0340 
0341 inline G4double G4Fragment::GetZ()  const
0342 {
0343   return static_cast<G4double>(theZ);
0344 }
0345 
0346 inline G4double G4Fragment::GetA() const
0347 {
0348   return static_cast<G4double>(theA);
0349 }
0350 
0351 inline void G4Fragment::SetZ(const G4double value)
0352 {
0353   theZ = G4lrint(value);
0354 }
0355 
0356 inline void G4Fragment::SetA(const G4double value)
0357 {
0358   theA = G4lrint(value);
0359 }
0360 
0361 inline G4int G4Fragment::GetNumberOfExcitons()  const
0362 {
0363   return numberOfParticles + numberOfHoles;
0364 }
0365 
0366 inline G4int G4Fragment::GetNumberOfParticles()  const
0367 {
0368   return numberOfParticles;
0369 }
0370 
0371 inline G4int G4Fragment::GetNumberOfCharged()  const
0372 {
0373   return numberOfCharged;
0374 }
0375 
0376 inline 
0377 void G4Fragment::SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
0378 {
0379   numberOfParticles = valueTot;
0380   numberOfCharged = valueP;
0381   if(valueTot < valueP)  { 
0382     NumberOfExitationWarning("SetNumberOfExcitedParticle"); 
0383   }
0384 }
0385 
0386 inline G4int G4Fragment::GetNumberOfHoles()  const
0387 {
0388   return numberOfHoles;
0389 }
0390 
0391 inline G4int G4Fragment::GetNumberOfChargedHoles()  const
0392 {
0393   return numberOfChargedHoles;
0394 }
0395 
0396 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
0397 {
0398   numberOfHoles = valueTot;
0399   numberOfChargedHoles = valueP;
0400   if(valueTot < valueP)  { 
0401     NumberOfExitationWarning("SetNumberOfHoles"); 
0402   }
0403 }
0404 
0405 inline void G4Fragment::SetNumberOfParticles(G4int value)
0406 {
0407   numberOfParticles = value;
0408 }
0409 
0410 inline void G4Fragment::SetNumberOfCharged(G4int value)
0411 {
0412   numberOfCharged = value;
0413   if(value > numberOfParticles)  { 
0414     NumberOfExitationWarning("SetNumberOfCharged"); 
0415   }
0416 }
0417 
0418 inline G4int G4Fragment::GetNumberOfElectrons() const
0419 {
0420   return numberOfShellElectrons;
0421 }
0422 
0423 inline void G4Fragment::SetNumberOfElectrons(G4int value)
0424 {
0425   numberOfShellElectrons = value;
0426 }
0427 
0428 inline G4int G4Fragment::GetCreatorModelID() const
0429 {
0430   return creatorModel;
0431 }
0432 
0433 inline void G4Fragment::SetCreatorModelID(G4int value)
0434 {
0435   creatorModel = value;
0436 }
0437 
0438 inline G4double G4Fragment::GetSpin() const
0439 {
0440   return spin;
0441 }
0442 
0443 inline void G4Fragment::SetSpin(G4double value)
0444 {
0445   spin = value;
0446 }
0447 
0448 inline G4bool G4Fragment::IsLongLived() const
0449 {
0450   return isLongLived;
0451 }
0452 
0453 inline void G4Fragment::SetLongLived(G4bool value)
0454 {
0455   isLongLived = value;
0456 }
0457 
0458 inline G4int G4Fragment::GetFloatingLevelNumber() const
0459 {
0460   return xLevel;
0461 }
0462 
0463 inline void G4Fragment::SetFloatingLevelNumber(G4int value)
0464 {
0465   xLevel = value;
0466 }
0467 
0468 inline 
0469 const G4ParticleDefinition* G4Fragment::GetParticleDefinition(void) const
0470 {
0471   return theParticleDefinition;
0472 }
0473 
0474 inline void G4Fragment::SetParticleDefinition(const G4ParticleDefinition * p)
0475 {
0476   theParticleDefinition = p;
0477 }
0478 
0479 inline G4double G4Fragment::GetCreationTime() const 
0480 {
0481   return theCreationTime;
0482 }
0483 
0484 inline void G4Fragment::SetCreationTime(G4double time)
0485 {
0486   theCreationTime = time;
0487 }
0488 
0489 inline G4NuclearPolarization* G4Fragment::NuclearPolarization()
0490 {
0491   return thePolarization;
0492 }
0493 
0494 inline G4NuclearPolarization* G4Fragment::GetNuclearPolarization() const
0495 {
0496   return thePolarization;
0497 }
0498 
0499 inline void G4Fragment::SetNuclearPolarization(G4NuclearPolarization* p)
0500 {
0501   thePolarization = p;
0502 }
0503 
0504 #endif