G4CascadeCoalescence.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // G4CascadeCoalescence:  Factory model for final-state interactions to
00027 //   produce light ions from cascade nucleons.  The algorithm implemented
00028 //   here is descirbed in Section 2.3 of the LAQGSM documentation (p. 11-12)
00029 //   [http://lib-www.lanl.gov/la-pubs/00818645.pdf].
00030 //
00031 // 20110917  Michael Kelsey
00032 // 20110920  M. Kelsey -- Use environment variables to set momentum cuts for tuning,
00033 //           replace polymorphic argument lists with use of "ClusterCandidate"
00034 
00035 #ifndef G4CASCADE_COALESCENCE_HH
00036 #define G4CASCADE_COALESCENCE_HH
00037 
00038 #include "globals.hh"
00039 #include "G4InuclNuclei.hh"
00040 #include "G4LorentzVector.hh"
00041 #include <vector>
00042 #include <set>
00043 
00044 class G4CollisionOutput;
00045 class G4InuclElementaryParticle;
00046 
00047 
00048 class G4CascadeCoalescence {
00049 public:
00050   G4CascadeCoalescence(G4int verbose=0);
00051   virtual ~G4CascadeCoalescence();
00052 
00053   // Final state particle list is modified directly
00054   void FindClusters(G4CollisionOutput& finalState);
00055 
00056   void setVerboseLevel(G4int verbose) { verboseLevel = verbose; }
00057 
00058 private:
00059   typedef std::vector<size_t> ClusterCandidate; // Indices of constituents
00060 
00061   G4int verboseLevel;                           // Control diagnostic messages
00062 
00063   static const G4double dpMaxDoublet;           // Relative momenta for clusters
00064   static const G4double dpMaxTriplet;
00065   static const G4double dpMaxAlpha;
00066 
00067   std::vector<ClusterCandidate> allClusters;    // List of candidates found
00068   std::set<size_t> triedClusters;               // Hashes of combinatorics
00069   std::set<size_t> usedNucleons;                // List of converted nucleons
00070 
00071   G4CollisionOutput* thisFinalState;            // Pointers to current event
00072   const std::vector<G4InuclElementaryParticle>* thisHadrons;
00073 
00074   ClusterCandidate thisCluster;                 // Reusable buffer for attempts
00075   G4InuclNuclei thisLightIon;                   // Reusable construction buffer
00076 
00077   // Processing stages -- search, construct, cleanup
00078   void selectCandidates();
00079   void createNuclei();
00080   void removeNucleons();
00081 
00082   // Do combinatorics of given nucleons to make candidates
00083   void tryClusters(size_t idx1, size_t idx2);
00084   void tryClusters(size_t idx1, size_t idx2, size_t idx3);
00085   void tryClusters(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
00086 
00087   // Create cluster candidate with listed indices
00088   void fillCluster(size_t idx1, size_t idx2);
00089   void fillCluster(size_t idx1, size_t idx2, size_t idx3);
00090   void fillCluster(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
00091 
00092   // Convert cluster to hash index (for combinatoric reduction)
00093   size_t clusterHash(const ClusterCandidate& clus) const;
00094 
00095   // Check if candidate cluster has already been evaluated
00096   bool clusterTried(const ClusterCandidate& clus) const {
00097     return triedClusters.find(clusterHash(clus)) != triedClusters.end();
00098   }
00099 
00100   // Check if indexed nucleon is already in a cluster
00101   bool nucleonUsed(size_t idx) const {
00102     return usedNucleons.find(idx) != usedNucleons.end();
00103   }
00104 
00105   // Evaluate conditions for cluster to form light ion
00106   bool allNucleons(const ClusterCandidate& clus) const;
00107   bool goodCluster(const ClusterCandidate& clus) const;
00108   G4int clusterType(const ClusterCandidate& aCluster) const;
00109 
00110   // Extract hadron from final state list
00111   const G4InuclElementaryParticle& getHadron(size_t idx) const {
00112     return (*thisHadrons)[idx];
00113   }
00114 
00115   // Convert candidate nucleon set into output nucleus (true == success)
00116   bool makeLightIon(const ClusterCandidate& aCluster);
00117 
00118   // Kinematics for cluster evaluations
00119   G4LorentzVector getClusterMomentum(const ClusterCandidate& aCluster) const;
00120 
00121   G4double maxDeltaP(const ClusterCandidate& aCluster) const;
00122 
00123   // Report cluster arguments for validation
00124   void reportArgs(const G4String& name, const ClusterCandidate& clus) const;
00125   void reportResult(const G4String& name, const G4InuclNuclei& nucl) const;
00126 };
00127 
00128 #endif  /* G4CASCADE_COALESCENCE_HH */

Generated on Mon May 27 17:47:49 2013 for Geant4 by  doxygen 1.4.7