Back to home page

EIC code displayed by LXR

 
 

    


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

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 // 20100315  M. Kelsey -- Remove "using" directory and unnecessary #includes.
0028 // 20100413  M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
0029 // 20100517  M. Kelsey -- Inherit from common base class, make other colliders
0030 //      simple data members
0031 // 20100617  M. Kelsey -- Make G4NucleiModel a data member, instead of
0032 //      creating and deleting on every cycle.
0033 // 20100623  M. Kelsey -- Undo change from 0617.  G4NucleiModel not reusable.
0034 // 20100714  M. Kelsey -- Switch to new G4CascadeColliderBase class
0035 // 20100716  M. Kelsey -- Eliminate inter_case; use base-class functionality,
0036 //      add function to compute recoil nuclear mass on the fly
0037 // 20100720  M. Kelsey -- Make EPCollider pointer member
0038 // 20100722  M. Kelsey -- Move cascade output buffers to .hh file
0039 // 20100728  M. Kelsey -- Move G4NucleiModel here, as pointer member
0040 // 20100907  M. Kelsey -- Add new "makeResidualFragment" to create
0041 //      G4InuclNuclei at current stage of cascade
0042 // 20100909  M. Kelsey -- Drop makeResidualFragment(), getResidualMass() and
0043 //      local G4InuclNuclei object, replace with new RecoilMaker.
0044 //      Move goodCase() to RecoilMaker.
0045 // 20100916  M. Kelsey -- Add functions to handle trapped particles, and to
0046 //      decay hyperons.
0047 // 20110224  M. Kelsey -- Add ::rescatter() function which takes a list of
0048 //      pre-existing secondaries as input.  Split ::collide() into
0049 //      separate utility functions.  Move cascade parameters to static
0050 //      data members.  Add setVerboseLevel().
0051 // 20110303  M. Kelsey -- Add more cascade functions to support rescattering
0052 // 20110304  M. Kelsey -- Modify rescatter to use original Propagate() input
0053 // 20110316  M. Kelsey -- Add function to do G4KineticTrack conversion, decay
0054 //      rescattering resonances in situ.
0055 // 20110324  M. Kelsey -- Add list of nucleon hit locations for rescatter().
0056 // 20110721  M. Kelsey -- Drop decayTrappedParticle(G4KineticTrack*).
0057 // 20110722  M. Kelsey -- Deprecate "output_particles" list in favor of using
0058 //      output directly (will help with pre-cascade issues).
0059 // 20110729  M. Kelsey -- Replace convertKineticToCascade() to reduce churn.
0060 // 20110801  M. Kelsey -- Add local target buffers for rescattering, to avoid
0061 //      memory leak.
0062 // 20110919  M. Kelsey -- Add optional final-state clustering
0063 // 20130304  M. Kelsey -- Add new G4CascadeHistory for cacasde structure reporting
0064 // 20130620  Address Coverity complaint about missing copy actions
0065 // 20141204  M. Kelsey -- Add function to test for non-interacting particles
0066 
0067 #ifndef G4INTRA_NUCLEI_CASCADER_HH
0068 #define G4INTRA_NUCLEI_CASCADER_HH
0069 
0070 #include "G4CascadeColliderBase.hh"
0071 #include "G4CollisionOutput.hh"
0072 #include "G4ThreeVector.hh"
0073 #include <vector>
0074 
0075 class G4CascadParticle;
0076 class G4CascadeCoalescence;
0077 class G4CascadeHistory;
0078 class G4CascadeRecoilMaker;
0079 class G4ElementaryParticleCollider;
0080 class G4InuclElementaryParticle;
0081 class G4InuclParticle;
0082 class G4KineticTrack;
0083 class G4KineticTrackVector;
0084 class G4NucleiModel;
0085 class G4V3DNucleus;
0086 
0087 
0088 class G4IntraNucleiCascader : public G4CascadeColliderBase {
0089 public:
0090   G4IntraNucleiCascader();
0091   virtual ~G4IntraNucleiCascader();
0092 
0093   void collide(G4InuclParticle* bullet, G4InuclParticle* target,
0094            G4CollisionOutput& globalOutput);
0095 
0096   // For use with Propagate to preload a set of secondaries
0097   void rescatter(G4InuclParticle* bullet, G4KineticTrackVector* theSecondaries,
0098          G4V3DNucleus* theNucleus, G4CollisionOutput& globalOutput);
0099 
0100   void setVerboseLevel(G4int verbose=0);
0101 
0102 private:
0103   static const G4int itry_max;      // Maximum number of attempts
0104   static const G4int reflection_cut;    // Maximum number of reflections
0105   static const G4double small_ekin; // Tolerance for round-off zero
0106   static const G4double quasielast_cut; // To recover elastic scatters
0107 
0108 protected:
0109   G4bool initialize(G4InuclParticle* bullet, G4InuclParticle* target);
0110 
0111   void newCascade(G4int itry);      // Clear buffers for next attempt
0112   void setupCascade();          // Fill cascade using nuclear model
0113   void generateCascade();       // Track secondaries through nucleus
0114   G4bool finishCascade();       // Clean up output, check consistency
0115 
0116   void finalize(G4int itry,         // Transfer final state for return
0117         G4InuclParticle* bullet, G4InuclParticle* target,
0118         G4CollisionOutput& globalOutput);
0119 
0120   G4InuclParticle* createTarget(G4V3DNucleus* theNucleus);
0121 
0122   // Functions to transfer input high-energy cascade for propagation
0123   void preloadCascade(G4V3DNucleus* theNucleus,
0124               G4KineticTrackVector* theSecondaries);
0125   void copyWoundedNucleus(G4V3DNucleus* theNucleus);
0126   void copySecondaries(G4KineticTrackVector* theSecondaries);
0127   void processSecondary(const G4KineticTrack* aSecondary);
0128   void releaseSecondary(const G4KineticTrack* aSecondary);
0129 
0130   // Functions to handle, e.g., low-energy hyperons stuck inside potential
0131   void processTrappedParticle(const G4CascadParticle& trapped);
0132   void decayTrappedParticle(const G4CascadParticle& trapped);
0133 
0134   // Test if particle is able to interact in nucleus
0135   G4bool particleCanInteract(const G4CascadParticle& cpart) const;
0136 
0137 private: 
0138   G4NucleiModel* model;
0139   G4ElementaryParticleCollider* theElementaryParticleCollider;
0140   G4CascadeRecoilMaker* theRecoilMaker;
0141   G4CascadeCoalescence* theClusterMaker;
0142   G4CascadeHistory* theCascadeHistory;
0143 
0144   // Buffers and parameters for cascade attempts
0145   G4InuclNuclei* tnuclei;       // Target nucleus (must be non-zero)
0146   G4InuclNuclei* bnuclei;       // Non-zero if ion-ion collision
0147   G4InuclElementaryParticle* bparticle; // Non-zero if hadron-ion collision
0148 
0149   G4double minimum_recoil_A;        // Require fragment with this mass
0150   G4double coulombBarrier;
0151 
0152   // Buffers for creation (and reuse) of rescattering targets
0153   G4InuclNuclei* nucleusTarget;
0154   G4InuclElementaryParticle* protonTarget;
0155 
0156   // Buffers for collecting result of cascade (reset on each iteration)
0157   G4CollisionOutput output;
0158   std::vector<G4CascadParticle> cascad_particles;
0159   std::vector<G4CascadParticle> new_cascad_particles;
0160   G4ExitonConfiguration theExitonConfiguration;
0161 
0162   std::vector<G4ThreeVector> hitNucleons;   // Nucleons hit before rescatter
0163 
0164 private:
0165   // Copying of modules is forbidden
0166   G4IntraNucleiCascader(const G4IntraNucleiCascader&);
0167   G4IntraNucleiCascader& operator=(const G4IntraNucleiCascader&);
0168 };        
0169 
0170 #endif /* G4INTRA_NUCLEI_CASCADER_HH */