TVenus


class description - source file - inheritance tree

class TVenus : public TGenerator


    public:
TVenus TVenus() TVenus TVenus(char* choice, Int_t numberOfEvents, Int_t zprojectile, Int_t aprojectile, Int_t ztarget, Int_t atarget, Double_t incidentMomentum, Bool_t labSys, Bool_t lastGeneration, Double_t impactParameterMin, Double_t impactParameterMax, Double_t phiMin, Double_t phiMax) TVenus TVenus(const TVenus&) virtual void ~TVenus() static TClass* Class() virtual void GenerateEvent(Option_t* option) virtual Int_t ImportParticles(TClonesArray* particles, Option_t* option) virtual TObjArray* ImportParticles(Option_t* option) virtual TClass* IsA() const virtual void SetAntiQuarkColourExchange(Int_t iaqu = 1) virtual void SetAveragep_T(Float_t ptf = 0.45) virtual void SetBaryonRadius(Float_t radbar = 0.63) virtual void SetBergerJaffeTheta(Float_t themas = 0.51225) virtual void SetBreakingProcedureOption(Int_t iopbrk = 1) virtual void SetCascadeStable(Bool_t stable) virtual void SetCentralPointInEvolution(Int_t iwcent = 0) virtual void SetClean(Int_t nclean = 0) virtual void SetCMToLabTransformation(Int_t labsys = 1) virtual void SetCollisionTrigger(Int_t ko1ko2 = 9999) virtual void SetCutOffForKmaxor(Int_t kutdiq = 4) virtual void SetDecayTime(Float_t taunll = 1.0) virtual void SetDensityDistributionRange(Float_t fctrmx = 10.0) virtual void SetDiffractiveValenceQuarkFrac(Float_t pvalen = 0.30) virtual void SetDualPartonModel(Int_t idpm = 0) virtual void SetEntropyCalculated(Int_t ientro = 2) virtual void SetEntropyOption(Int_t iopent = 5, Float_t uentro = 4., Int_t kentro = 10000) virtual void SetEventPrint(Int_t ishevt = 0) virtual void SetGaussDistributionRange(Float_t gaumx = 8.0) virtual void SetGribovReggeCrossSecWeight(Float_t gricel = 1.5) virtual void SetGribovReggeDelta(Float_t gridel = 0.07) virtual void SetGribovReggeGamma(Float_t grigam = 3.64*0.04) virtual void SetGribovReggeRSquared(Float_t grirsq = 3.56*0.04) virtual void SetGribovReggeSlope(Float_t grislo = 0.25*0.04) virtual void SetHardCoreDistance(Float_t core = 0.8) virtual void SetHeavyFlavoursSpinProb(Float_t pspinh = 0.75) virtual void SetInelasticProtonProtonCrossSec(Float_t sigppi = -1.0) virtual void SetInitialRandomSeed(Double_t seedi = 0.0) virtual void SetInteractionMass(Float_t amsiac = 0.8) virtual void SetIsoSpinProb(Float_t pispn = 0.50) virtual void SetJ_PsiEvolutionTime(Float_t taumx = 20.0) virtual void SetJ_PsiEvolutionTimeSteps(Int_t nsttau = 100) virtual void SetJ_PsiNucleonCrossSec(Float_t sigj = 0.2) virtual void SetJFRADESuppression(Int_t ifrade = 1) virtual void SetJIntaOption(Int_t iojint = 2) virtual void SetKShortKLongStable(Bool_t stable) virtual void SetLambdaStable(Bool_t stable) virtual void SetLightFlavoursSpinProb(Float_t pspinl = 0.50) virtual void SetMaxImpact(Float_t bmaxim = 10000.) virtual void SetMaxNumberOfCollisions(Int_t ncolmx = 10000) virtual void SetMaxNumberOfValenceQuarks(Int_t neqmx = 5) virtual void SetMaxResonanceSpin(Int_t maxres = 99999) virtual void SetMeanTransverseQuarkMomentum(Float_t ptq1 = 0.260, Float_t ptq2 = 0.0, Float_t ptq3 = 0.0) virtual void SetMesonRadius(Float_t radmes = 0.40) virtual void SetMinImpact(Float_t bminim = 0.) virtual void SetMinimumEnergyOption(Int_t iopenu = 1) virtual void SetMinNumberOfValenceQuarks(Int_t neqmn = -5) virtual void SetMinTimeInEvolution(Float_t wtmini = -3.0) virtual void SetMomentumRescaling(Int_t irescl = 1) virtual void SetMuonAngle(Float_t angmue = 3.9645/360.0*2*3.14159) virtual void SetMuonEnergy(Float_t elepto = 26.24) virtual void SetNueEnergy(Float_t elepti = 43.00) virtual void SetOmegaStable(Bool_t stable) virtual void SetOscillatorQuantum(Float_t omega = .500) virtual void Setp_T_Distribution(Int_t ioptf = 1) virtual void Setp_TDistributionRange(Float_t ptmx = 6.0) virtual void SetPhaseSpace(Float_t delmss = 0.300) virtual void SetPiZeroStable(Bool_t stable) virtual void SetPrintMarks(Int_t ipagi = 0) virtual void SetPrintOption(Int_t ish = 0) virtual void SetPrintOptionAmprif(Float_t amprif = 0.0) virtual void SetPrintOptionDeleps(Float_t deleps = 1.0) virtual void SetPrintOptionDelvol(Float_t delvol = 1.0) virtual void SetPrintSubOption(Int_t ishsub = 0) virtual void SetProjectileDiffractiveProb(Float_t wproj = 0.32) virtual void SetQQ_QQbarProbability(Float_t pdiqua = 0.08) virtual void SetQuarkp_TDistributionOption(Int_t ioptq = 2) virtual void SetRapidityUpperLimit(Float_t ymximi = 2.0) virtual void SetReactionTime(Float_t taurea = 1.5) virtual void SetResonanceStable(Bool_t stable = kFALSE) virtual void SetRhoPhiRatio(Float_t rhophi = 0.5) virtual void SetSeaProbability(Float_t prosea = -1.0) virtual void SetSeaRatio(Float_t rstras = 0.0) virtual void SetSemihardCutOff(Float_t pth = 1.0) virtual void SetSemihardInteractionProb(Float_t phard = -1.0) virtual void SetSigmaStable(Bool_t stable) virtual void SetsMass(Float_t smas = 0.0) virtual void SetSpaceTimeEvolution(Int_t ispall = 1) virtual void SetSpaceTimeEvolutionMinTau(Float_t taumin = 1.) virtual void SetssMass(Float_t ssmas = 0.0) virtual void SetStorage(Int_t istore) virtual void SetStoreOnlyStable(Int_t istmax = 0) virtual void SetStringDecayParameter(Float_t parea = 0.20) virtual void SetStringTension(Float_t tensn = 1.0) virtual void SetStructureFunctionCutOffMass(Float_t cutmss = 0.001) virtual void SetStructureFunctionSeaValence(Float_t cutmsq = 2.0) virtual void SetTargetDiffractiveProb(Float_t wtarg = 0.32) virtual void SetTauSteps(Int_t numtau = 86) virtual void SetThresholdResonanceToString(Float_t delrem = 1.0) virtual void SetTimeStepInEvolution(Float_t wtstep = 1.0) virtual void SetTryAgain(Int_t ntrymx = 10) virtual void SetU_D_QuarkProductionProb(Float_t pud = 0.455) virtual void SetusMass(Float_t usmas = 0.0) virtual void SetuuMass(Float_t uumas = 0.0) virtual void SetVersionNumber(Int_t iversn = 521) virtual void Show() virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
Bool_t fAllStable Int_t fAProjectile Double_t fAveragePt Int_t fATarget Double_t fBaryonRadius Int_t fBreaking Bool_t fCascadeStable Int_t fCentralPoint TString fChoice Int_t fClean Int_t fCollisionTrigger Int_t fColourExchange Double_t fCore Int_t fCutKmaxor Double_t fCutMass Double_t fCutSemi Double_t fDecayTime Double_t fDensityRange Int_t fDualParton Int_t fEntropy1 Double_t fEntropy2 Int_t fEntropy3 Int_t fEntropyCalc Double_t fGaussianRange Double_t fGribovCross Double_t fGribovDelta Double_t fGribovGamma Double_t fGribovR2 Double_t fGribovSlope Double_t fImpactParameterMin Double_t fImpactParameterMax Double_t fInitialSeed Double_t fInteractionMass Double_t fInterval Int_t fIstore Int_t fJfradeSup Int_t fJintalOpt Double_t fJ_PsiEvolution Int_t fJ_PsiSteps Bool_t fKStable Bool_t fLabSys Bool_t fLambdaStable Bool_t fLastGeneration Double_t fMassTheta Int_t fMaxCollisions Int_t fMaxSpin Int_t fMaxValence Double_t fMeanQMomentum1 Double_t fMeanQMomentum2 Double_t fMeanQMomentum3 Double_t fMesonRadius Int_t fMinEnergy Double_t fMinTime Int_t fMinValence Double_t fMuonAngle Double_t fMuonEnergy Double_t fNueEnergy Int_t fNumberOfEvents Double_t fIncidentMomentum Bool_t fOmegaStable Double_t fOsciQuantum Double_t fPIsospin Bool_t fPiZeroStable Double_t fPhaseSpace Double_t fPhiMin Double_t fPhiMax Double_t fPDiQuark Double_t fPPInelastic Int_t fPrintEvent Int_t fPrintMarks Int_t fPrintOption Int_t fPrintOption1 Int_t fPrintOption2 Int_t fPrintOption3 Int_t fPrintSubOption Double_t fProjDiffProb Double_t fPSpinLight Double_t fPSpinHeavy Int_t fPTDistribution Double_t fPtRange Double_t fPud Int_t fQuarkPt Double_t fReacTime Int_t fRescale Double_t fResThreshold Int_t fRetries Double_t fRhoPhi Double_t fSeaProbability Double_t fSeaRatio Double_t fSeaValenceCut Double_t fSemihard Double_t fSigmaJ_Psi Bool_t fSigmaStable Int_t fSpaceTime Double_t fSpaceTimeStep Double_t fssQuarkMass Double_t fSQuarkMass Double_t fStringDecay Double_t fStringTension Double_t fTargDiffProb Double_t fTauMin Double_t fTauMax Int_t fTauSteps Bool_t fUpdate Double_t fuuQuarkMass Double_t fusQuarkMass Double_t fValenceFrac Int_t fVersion Int_t fZProjectile Int_t fZTarget

Class Description

                                                                          
 TVenus                                                                   
                                                                          
 TVenus is an interface class between the event generator VENUS and       
 the ROOT system. The current implementation is based on VENUS 5.21       
                                                                          
 Authors of VENUS:                                                        
           Klaus WERNER                                                   
           Eugen FURLER (upgrade 5.04)                                    
           Michael HLADIK (upgrades 5.05, 5.13, 5.15, 5.17)               
                                                                          
  Laboratoire Subatech, Universite de Nantes - IN2P3/CNRS - EMN           
  4 rue Alfred Kastler, 44070 Nantes Cedex 03, France                     
                                                                          
  Email: <last name>@nanhp3.in2p3.fr                                      
                                                                          
  VENUS is a Monte Carlo procedure to simulate hadronic interactions at   
  ultrarelativistic energies (hadron-hadron, hadron-nucleus, nucleus-     
  nucleus scattering).                                                    
                                                                          
  VENUS is based on Gribov-Regge theory (of multiple Pomeron exchange)    
  and classical relativistic string dynamics. A detailed description can  //  
  be found in a review article, published in Physics Reports 232 (1993)   
  pp. 87-299.                                                             
                                                                          
                                                                          


TVenus():TGenerator("Venus","Venus")
  Event generator Venus default constructor


TVenus(char *choice , Int_t /*numberOfEvents*/, Int_t zprojectile , Int_t aprojectile, Int_t ztarget , Int_t atarget, Double_t incidentMomentum, Bool_t labSys , Bool_t lastGeneration, Double_t impactParameterMin, Double_t impactParameterMax, Double_t phiMin , Double_t phiMax) :TGenerator("Venus","Venus")
  Event generator Venus normal constructor


~TVenus()
  Event generator VENUS default destructor


void GenerateEvent(Option_t *option)
  Generate next event


TObjArray* ImportParticles(Option_t *option)
  Overloaded primary creation method. The event generator does
  not use the HEPEVT common block, and has not the PDG numbering scheme.


Int_t ImportParticles(TClonesArray *particles, Option_t *option)
  Default primary creation method. It reads the /HEPEVT/ common block which
  has been filled by the GenerateEvent method. If the event generator does
  not use the HEPEVT common block, This routine has to be overloaded by
  the subclasses.
  The function loops on the generated particles and store them in
  the TClonesArray pointed by the argument particles.
  The default action is to store only the stable particles (ISTHEP = 1)
  This can be demanded explicitly by setting the option = "Final"
  If the option = "All", all the particles are stored.


void SetVersionNumber(Int_t iversn)
  Set VENUS version number, the default is version 5.21


void SetU_D_QuarkProductionProb(Float_t pud)
  Set VENUS ud quark production probability


void SetQQ_QQbarProbability(Float_t pdiqua)
  Set VENUS qq - qqbar probability


void SetLightFlavoursSpinProb(Float_t pspinl)
  Set VENUS light flavour spin probability


void SetHeavyFlavoursSpinProb(Float_t pspinh)
  Set VENUS heavy flavour spin probability


void SetIsoSpinProb(Float_t pispn)
  Set VENUS isospin probability


void Setp_T_Distribution(Int_t ioptf)
  Set VENUS transverse momentum distribution
  ioptf = 1 -> exponential
  ioptf = 2 -> gaussian


void SetAveragep_T(Float_t ptf)
  Set VENUS average transverse momentum


void SetStringTension(Float_t tensn)
  Set VENUS string tension


void SetStringDecayParameter(Float_t parea)
  Set VENUS string decay parameter


void SetThresholdResonanceToString(Float_t delrem)
  Set VENUS threshold for resonance to string


void SetCutOffForKmaxor(Int_t kutdiq)
  Set VENUS cut off for Kmaror


void SetBreakingProcedureOption(Int_t iopbrk)
  Set VENUS breaking procedure option
  iopbrk = 1 -> amor
  iopbrk = 2 -> samba


void SetQuarkp_TDistributionOption(Int_t ioptq)
  Set VENUS quark p_T distribution option
  ioptq = 1 -> exponential
  ioptq = 2 -> gaussian
  ioptq = 3 -> power


void SetMeanTransverseQuarkMomentum(Float_t ptq1, Float_t ptq2, Float_t ptq3)
  Set VENUS mean quark transverse momentum
  with the energy dependence of mean transverse momentum of quarks:
  ptq1+ptq2*alog(e)+ptq3*alog(e)**2
  with e=sqrt(s)


void SetSemihardInteractionProb(Float_t phard)
  Set VENUS semihard interaction probability (not used if negative)


void SetSemihardCutOff(Float_t pth)
  Set VENUS cutoff parameter for p_t distr for semihard interactions


void SetSeaRatio(Float_t rstras)
  Set VENUS ratio of strang sea quark contents over up sea contents


void SetProjectileDiffractiveProb(Float_t wproj)
  Set VENUS projectile diffractive probability


void SetTargetDiffractiveProb(Float_t wtarg)
  Set VENUS target diffractive probability


void SetStructureFunctionSeaValence(Float_t cutmsq)
  Set VENUS effective cutoff mass in structure fcts for sea-valence ratio


void SetStructureFunctionCutOffMass(Float_t cutmss)
  Set VENUS effective cutoff mass in structure fcts for sea-valence ratio


void SetDiffractiveValenceQuarkFrac(Float_t pvalen)
  Set VENUS valence quark fraction in case of diffractive interaction


void SetPhaseSpace(Float_t delmss)
  Set VENUS phase space parameter


void SetGribovReggeGamma(Float_t grigam)
  Set VENUS gribov-regge-theory gamma


void SetGribovReggeRSquared(Float_t grirsq)
  Set VENUS gribov-regge-theory r^2


void SetGribovReggeDelta(Float_t gridel)
  Set VENUS gribov-regge-theory delta


void SetGribovReggeSlope(Float_t grislo)
  Set VENUS gribov-regge-theory slope


void SetGribovReggeCrossSecWeight(Float_t gricel)
  Set VENUS gribov-regge-theory cross section weight


void SetHardCoreDistance(Float_t core)
  Set VENUS hard core distance


void SetJ_PsiNucleonCrossSec(Float_t sigj)
  Set VENUS J/Psi nucleon cross section


void SetReactionTime(Float_t taurea)
  Set VENUS reaction time


void SetBaryonRadius(Float_t radbar)
  Set VENUS baryon radius


void SetMesonRadius(Float_t radmes)
  Set VENUS meson radius


void SetInteractionMass(Float_t amsiac)
  Set VENUS interaction mass


void SetJIntaOption(Int_t iojint)
  Set VENUS option to call jinta1 (1) or jinta2 (2)


void SetPrintOptionAmprif(Float_t amprif)
  Set VENUS print option


void SetPrintOptionDelvol(Float_t delvol)
  Set VENUS print option


void SetPrintOptionDeleps(Float_t deleps)
  Set VENUS print option


void SetEntropyOption(Int_t iopent, Float_t uentro, Int_t kentro)
  Set VENUS options for entropy calculation
  option for entropy calculation
  iopent = 0 -> zero entropy
  iopent = 1 -> oscillator model (0 for k.le.uentro)
  iopent = 2 -> fermi gas w const volume (0 for k.le.uentro)
  iopent = 3 -> fermi gas w const density (0 for k.le.uentro)
  iopent = 4 -> fermi gas w const vol - new (0 for k.le.uentro)
  iopent = 5 -> resonance gas (hagedorn) (0 for u.le.uentro)


void SetDecayTime(Float_t taunll)
  Set VENUS decay time of comoving frame


void SetOscillatorQuantum(Float_t omega)
  Set VENUS oscillator quantum


void SetSpaceTimeEvolutionMinTau(Float_t taumin)
  Set VENUS minimum tau for space-time evolution


void SetTauSteps(Int_t numtau)
  Set VENUS number of tau steps


void Setp_TDistributionRange(Float_t ptmx)
  Set VENUS transverse momentum distribution range


void SetGaussDistributionRange(Float_t gaumx)
  Set VENUS gaussian distribution range


void SetDensityDistributionRange(Float_t fctrmx)
  Set VENUS density distribution range


void SetTryAgain(Int_t ntrymx)
  Set VENUS number of retries


void SetJ_PsiEvolutionTime(Float_t taumx)
  Set VENUS J/Psi evolution time


void SetJ_PsiEvolutionTimeSteps(Int_t nsttau)
  Set VENUS J/Psi evolution time steps


void SetMinimumEnergyOption(Int_t iopenu)
  Set VENUS minimum energy option
  iopent = 1 -> sum of hadron masses
  iopent = 2 -> bag model curve with minimum at nonzero strangen


void SetBergerJaffeTheta(Float_t themas)
  Set VENUS parameter theta in berger/jaffe mass formula


void SetSeaProbability(Float_t prosea)
  Set VENUS sea prob (if < 0 then calculated from structure fncts)


void SetInelasticProtonProtonCrossSec(Float_t sigppi)
  Set VENUS inelastic pp cross section [fm**2]
  if negative: calculated from gribov-regge-theory


void SetEntropyCalculated(Int_t ientro)
  Set VENUS entro() calculated (1) or from data (2)


void SetDualPartonModel(Int_t idpm)
  Set VENUS dual parton model (1) or not (else)


void SetAntiQuarkColourExchange(Int_t iaqu)
  Set VENUS antiquark color exchange (1) or not (0)


void SetMinNumberOfValenceQuarks(Int_t neqmn)
  Set VENUS minimum number of valence quarks


void SetMaxNumberOfValenceQuarks(Int_t neqmx)
  Set VENUS maximum number of valence quarks


void SetRapidityUpperLimit(Float_t ymximi)
  Set VENUS upper limit for rapidity interval for intermittency analysis


void SetClean(Int_t nclean)
  Set VENUS clean /cptl/ if nclean > 0 (every nclean_th time step)


void SetCMToLabTransformation(Int_t labsys)
  Set VENUS trafo from pp-cm into lab-system (1) or not (else)


void SetMaxNumberOfCollisions(Int_t ncolmx)
  Set VENUS maximum number of collisions


void SetMaxResonanceSpin(Int_t maxres)
  Set VENUS maximum resonance spin


void SetMomentumRescaling(Int_t irescl)
  Set VENUS momentum rescaling (1) or not (else)


void SetNueEnergy(Float_t elepti)
  Set VENUS nue energy


void SetMuonEnergy(Float_t elepto)
  Set VENUS muon energy


void SetMuonAngle(Float_t angmue)
  Set VENUS muon angle


void SetCollisionTrigger(Int_t ko1ko2)
  Set VENUS collision trigger (only coll between ko1 and ko2 are used)


void SetPrintOption(Int_t ish)
  Set VENUS print option
    ish=0      option for printout
        14 -> call uttima
        16 -> prints sea prob.
        18 -> sr jclude, no-phase-space droplets
        19 -> sr ainitl, call smassp
        21 -> creates histogram for sea distribution
        22 -> sr jfrade, msg after call utclea
        23 -> call jintfp
        24 -> call jintcl
        25 -> call jchprt
        26 -> call qgcpen (plot of energy and flavor density)
        27 -> call qnbcor
        29 -> call qgcpfl (plot of dispersion)
        90,91,92,93,94,95 -> more and more detailed messages.


void SetPrintSubOption(Int_t ishsub)
  Set VENUS print sub option
      ishsub=0   option for printout
      ishsub=ijmn, ij specifies location where ish=mn.
             ij=01 -> sr jcludc
             ij=02 -> sr jetgen
             ij=03 -> sr jfrade, starting before fragmentation
             ij=04 -> sr jdecay
             ij=05 -> sr jdecax
             ij=06 -> sr nucoll
             ij=07 -> sr nucoge+-
             ij=08 -> sr aastor
             ij=09 -> sr jfrade, starting after fragmentation
             ij=10 -> sr jfrade, starting before decay
             ij=11 -> sr jfrade, starting after interactions
             ij=12 -> sr jcentr, entro() in data format
             ij=13 -> sr jcentp
             ij=14 -> sr jdecax if droplet decay
             ij=15 -> sr jsplit
             ij=16 -> sr jfrade
             ij=17 -> sr racpro
             ij=18 -> sr utclea
             ij=19 -> sr jinta1, jinta2, after call utclea
             ij=20 -> sr jdecas
             ij=21 -> sr jdecas (without jdecax)
             ij=22 -> sr utcley
             ij=50 -> sr qgcnbi


void SetEventPrint(Int_t ishevt)
  Set VENUS ishevt != 0: for evt# != ishevt ish is set to 0


void SetPrintMarks(Int_t ipagi)
  Set VENUS print marks between whom ish is set to ish(init)


void SetMaxImpact(Float_t bmaxim)
  Set VENUS maximum beam impact parameter


void SetMinImpact(Float_t bminim)
  Set VENUS minimum beam impact parameter


void SetStoreOnlyStable(Int_t istmax)
  Set VENUS store only stable ptl (0) or also parents (1)


void SetInitialRandomSeed(Double_t seedi)
  Set VENUS initial random number seed


void SetJFRADESuppression(Int_t ifrade)
  Set VENUS suppression of calling jfrade (0). jfrade=fragm+decay+rescatt


void SetResonanceStable(Bool_t stable)
  Set VENUS all resonance decays suppressed


void SetKShortKLongStable(Bool_t stable)
  Set VENUS k_short/long (+-20) decays suppressed


void SetLambdaStable(Bool_t stable)
  Set VENUS lambda (+-2130) decays suppressed


void SetSigmaStable(Bool_t stable)
  Set VENUS sigma (+-1130,+-2230) decays suppressed


void SetCascadeStable(Bool_t stable)
  Set VENUS cascade (+-2330,+-1330) decays suppressed


void SetOmegaStable(Bool_t stable)
  Set VENUS omega (+-3331) decays suppressed


void SetPiZeroStable(Bool_t stable)
  Set VENUS pi0 (110) decays suppressed


void SetRhoPhiRatio(Float_t rhophi)
  Set VENUS rho/rho+phi ratio


void SetSpaceTimeEvolution(Int_t ispall)
  Set VENUS wspa: all ptls (1) or only interacting ptls (else)


void SetMinTimeInEvolution(Float_t wtmini)
  Set VENUS tmin in wspa


void SetTimeStepInEvolution(Float_t wtstep)
  Set VENUS t-step in wspa


void SetCentralPointInEvolution(Int_t iwcent)
  Set VENUS only central point (1) or longitudinal distr (else) in wspa


void SetsMass(Float_t smas)
  Set VENUS s quark mass


void SetuuMass(Float_t uumas)
  Set VENUS uu diquark mass


void SetusMass(Float_t usmas)
  Set VENUS us diquark mass


void SetssMass(Float_t ssmas)
  Set VENUS ss diquark mass


void SetStorage(Int_t istore)
  Set VENUS internal storage parameter


void Show()
  Shows the actual values of some common block variables as in
  the original VENUS command 'show'




Inline Functions


            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void Streamer(TBuffer& b)
               void StreamerNVirtual(TBuffer& b)
             TVenus TVenus(const TVenus&)


Author: Ola Nordmann 21/09/95
Last update: root/venus:$Name: $:$Id: TVenus.cxx,v 1.1.1.1 2000/05/16 17:00:48 rdm Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.