ATLFast


class description - source file - inheritance tree

class ATLFast : public TNamed


    public:
ATLFast ATLFast() ATLFast ATLFast(char* name, char* title = The ATLAS fast MonteCarlo) ATLFast ATLFast(ATLFast&) void ~ATLFast() Bool_t Bfield() virtual void Browse(TBrowser* b) TClass* Class() virtual void Clear(Option_t* option) ATLFClusterMaker* ClusterMaker() ATLFVirtualDisplay* Display() virtual void Draw(Option_t* option) ATLFElectronMaker* ElectronMaker() Int_t Event() virtual void FillClone() void FillTree() virtual void Finish() virtual void GetEvent(Int_t event = 1) Int_t GetVersion() Int_t GetVersionDate() virtual void Init() void InitChain(TChain* chain) virtual TClass* IsA() virtual Bool_t IsFolder() ATLFJetMaker* JetMaker() Int_t Luminosity() virtual Int_t Make(Int_t i = 0) ATLFMaker* Maker(char* name) TList* Makers() void MakeTree(char* name = T, char* title = ATLFast tree) ATLFMCMaker* MCMaker() ATLFMiscMaker* MiscMaker() Int_t Mode() ATLFMuonMaker* MuonMaker() virtual void Paint(Option_t* option) ATLFPhotonMaker* PhotonMaker() virtual void PrintInfo() Int_t Run() void SetBfield(Bool_t field = 1) virtual void SetDefaultParameters() void SetDisplay(ATLFVirtualDisplay* display) virtual void SetEvent(Int_t event = 1) void SetLuminosity(Int_t lumi = 1) virtual void SetMode(Int_t mode = 0) virtual void SetRun(Int_t run = 1) void SetSmearing(Bool_t val = 1) void SetSmearMuonFun(Int_t val = 2) void SetSmearMuonOpt(Int_t val = 3) void SetSUSYcodeLSP(Int_t val = 66) void SetTrackFinding(Bool_t val = 0) virtual void ShowMembers(TMemberInspector& insp, char* parent) Bool_t Smearing() Int_t SmearMuonFun() Int_t SmearMuonOpt() void SortDown(Int_t n, Float_t* a, Int_t* index, Bool_t down = kTRUE) virtual void Streamer(TBuffer& b) Int_t SUSYcodeLSP() Bool_t TrackFinding() ATLFTrackMaker* TrackMaker() TTree* Tree() ATLFTriggerMaker* TriggerMaker()

Data Members

private:
Int_t m_Version ATLFast version number Int_t m_VersionDate ATLFast version date Int_t m_Run Run number Int_t m_Event Event number Int_t m_Mode Run mode TTree* m_Tree Pointer to the Root tree TList* m_Makers List of Makers ATLFMCMaker* m_MCMaker Pointer to MCMaker ATLFClusterMaker* m_ClusterMaker Pointer to ClusterMaker ATLFElectronMaker* m_ElectronMaker Pointer to ElectronMaker ATLFMuonMaker* m_MuonMaker Pointer to MuonMaker ATLFPhotonMaker* m_PhotonMaker Pointer to PhotonMaker ATLFJetMaker* m_JetMaker Pointer to JetMaker ATLFTrackMaker* m_TrackMaker Pointer to TrackMaker ATLFTriggerMaker* m_TriggerMaker Pointer to TriggerMaker ATLFMiscMaker* m_MiscMaker Pointer to MiscMaker Int_t m_Luminosity Luminosity option (low=1, high=2) Bool_t m_Bfield B-field (on=1,off=0) Bool_t m_Smearing Smearing on=1, off=0 Int_t m_SmearMuonOpt Smear for muons: muon=1, track=2,combin=3 Int_t m_SmearMuonFun Smear for muons: exact=1, approxim=2 Int_t m_SUSYcodeLSP Code for SUSY LSP particle Int_t m_TrackFinding Track/finding on=1,off=0 ATLFHistBrowser m_HistBrowser Object to Browse Maker histograms ATLFBigBang m_BigBang !Object to Browse generated particles ATLFVirtualDisplay* m_Display !Pointer to Event display object

Class Description

                                                                      
 ATLFast                                                              
                                                                      
 Main class to control the ATLFast program.                           
                                                                      
 This class is a complete redesign of the Fortran program ATLfast.    
 This class is a general framework for programs that needs to:        
    - Initialise some parameters                                      
    - Loop on events                                                  
    - Print results and save histograms, etc                          
                                                                      
 The event processor ATLFast::Make loops on a list of Makers          
 where each maker performs some task on the event data and generates  
 results.                                                             
 New Makers can be inserted by a user without modifying this class.   
 Note that the order in which the Makers are called is the order      
 of insertion in the list of Makers.                                  
 Each Maker is responsible for creating its branch of the Tree.       
 The following table shows the list of makers currently implemented   
 The default option to Save the Maker info in the Tree is mentioned.  
                                                                      
    Maker name        Save in Tree                                    
    ==========        ============                                    
    MCMaker             NO                                            
    ClusterMaker        YES                                           
    ElectronMaker       YES                                           
    MuonMaker           YES                                           
    PhotonMaker         YES                                           
    JetMaker            YES                                           
    TrackMaker          NO                                            
    TriggerMaker        YES                                           
    MiscMaker           YES                                           
                                                                      
 Makers must derive from the base class ATLFMaker.                    
 ATLFMaker provides a common interface to all Makers.                 
 Each Maker is responsible for defining its own parameters and        
 histograms.                                                          
 Each Maker has its own list of histograms.                           
 Each Maker has an associated companion class corresponding to the    
 type of physics object reconstructed by the Maker.                   
 For example, ATLFClusterMaker creates ATLFCluster objects.           
              ATLFTriggerMaker creates one single ATLFTrigger object. 
 The pointer supporting the created object(s) is defined in ATLFMaker 
   m_Fruits may point to a single object (eg. ATLFTrigger) or to a    
           TClonesArray of objects (eg. ATLFCluster).                 
                                                                      
 The function ATLFast::Maketree must be called after the creation     
 of the ATLFast object to create a Root Tree.                         
                                                                      
 An example of main program/macro to use ATLFast is given below:      
========================================================================
void umain(Int_t nevents=100)
{
   gROOT->Reset();
   gSystem->Load("libatlfast.so");  // dynamically link the compiled shared library

   // Open the root output file
   TFile file("atlfast.root","recreate","ATLFast root file",2);

   ATLFast atlfast("atlfast");     // create main object to run atlfast

   User user;           // create an object of the User class defined in user.C

   atlfast.Init();      // Initialise event (maker histograms,etc)
   atlfast.MakeTree();  // Create the Root tree

   gROOT->LoadMacro("user.C");  // compile/interpret user file

   for (Int_t i=0; i<nevents; i++) {
      if (i%100 == 0) printf("In loop:%dn",i);
      atlfast.Make(i);       // Generate and reconstruct event
      user.FillHistograms(); // User has possibility to decide if store event here!
      atlfast.FillTree();
      atlfast.Clear();       // Clear all event lists
   }
   atlfast.Finish();

   // save objects in Root file
   atlfast.Write();  //save main atlfast object (and run parameters)
}
========================================================================
                                                                      
 This example illustrates how to:                                     
    - Load a shared library                                           
    - Open a Root file                                                
    - Initialise ATLFast                                              
    - Load some user code (interpreted)                               
      This user code may redefine some Maker parameters               
    - Make a loop on events                                           
    - Save histograms and the main ATLFast object and its Makers      
                                                                      
========================================================================
  An example of a User class is given below:                          
========================================================================

#ifndef user_H
#define user_H


                                                                      
 User                                                                 
                                                                      
 Example of a user class to perform user specific tasks when running  
 the ATLfast program.                                                 
                                                                      
 This class illustrates:                                              
   - How to set run parameters                                        
   - How to create and fill histograms                                
                                                                      


class TH1F;
class ATLFast;
class ATLFClusterMaker;
class ATLFPhotonMaker;

class User {

private:
   TH1F             *m_hist1;       //pointer to histogram
   TH1F             *m_hist2;       //pointer to histogram
   TH1F             *m_hist3;       //pointer to histogram
public:
               User();
   void        FillHistograms();
   void        SetRunParameters();

#endif
};

_________________________________________________________________________
User::User()
{
   SetRunParameters();  //change default parameters

         Create a few histograms
   m_hist1 = new TH1F("hist1","Number of tracks per event",100,0,100);
   m_hist2 = new TH1F("hist2","Number of clusters",100,0,100);
   m_hist3 = new TH1F("hist3","Number of isolated muons",20,0,20);
}

_________________________________________________________________________
void User::FillHistograms()
{
//   m_hist1.Fill(event->GetNtracks());
//   m_hist2.Fill(event->GetNclusters));
//   m_hist3.Fill(event->GetNIsoMuons());
}

_________________________________________________________________________
void User::SetRunParameters()
{
  // change ATLfast default parameters

   gATLFast->SetSmearMuonOpt(0);
   gATLFast->ClusterMaker()->SetGranBarrelEta(0.12);
   gATLFast->PhotonMaker()->SetMinPT(6.);
   gATLFast->TriggerMaker()->SetMuoEtaCoverage(2.8);

}
======================end of User class=================================



ATLFast() : TNamed("atlfast","The ATLAS fast simulation")

ATLFast(const char *name, const char *title) : TNamed(name,title)

~ATLFast()
   m_Makers->Delete();
   delete m_Makers;

void Browse(TBrowser *b)

void Clear(Option_t *option)
    Reset lists of event objects

void Draw(Option_t *option)
    Insert current event in graphics pad list

void GetEvent(Int_t event)
    Read event from Tree

void Init()
    Initialise makers

void Paint(Option_t *option)
    Paint ATLFast objects

void PrintInfo()
     Gives information about versions etc.

void FillTree()
  Fill the ROOT tree, looping on all active branches

void InitChain(TChain *chain)
  Initialize branch addresses for all makers in a TChain

void MakeTree(const char* name, const char*title)
  Create a ROOT tree
  Loop on all makers to create the Root branch (if any)

void SetDefaultParameters()

Int_t Make(Int_t i)

void FillClone()
 Fill Makers fruits clones

void Finish()
    Terminate a run
   place to make operations on histograms, normalization,etc.

void SortDown(Int_t n1, Float_t *a, Int_t *index, Bool_t down)
  sort the n1 elements of array a.
  In output the array index contains the indices of the sorted array.
  if down is false sort in increasing order (default is decreasing order)
   This is a translation of the CERNLIB routine sortzv (M101)
   based on the quicksort algorithm

void Streamer(TBuffer &R__b)
 Stream an object of class ATLFast.



Inline Functions


                      Int_t GetVersion()
                      Int_t GetVersionDate()
        ATLFVirtualDisplay* Display()
                     Bool_t IsFolder()
                     TList* Makers()
                 ATLFMaker* Maker(char* name)
               ATLFMCMaker* MCMaker()
          ATLFClusterMaker* ClusterMaker()
         ATLFElectronMaker* ElectronMaker()
             ATLFMuonMaker* MuonMaker()
           ATLFPhotonMaker* PhotonMaker()
              ATLFJetMaker* JetMaker()
            ATLFTrackMaker* TrackMaker()
          ATLFTriggerMaker* TriggerMaker()
             ATLFMiscMaker* MiscMaker()
                     TTree* Tree()
                      Int_t Run()
                      Int_t Event()
                      Int_t Mode()
                      Int_t Luminosity()
                     Bool_t Bfield()
                     Bool_t Smearing()
                      Int_t SmearMuonOpt()
                      Int_t SmearMuonFun()
                      Int_t SUSYcodeLSP()
                     Bool_t TrackFinding()
                       void SetDisplay(ATLFVirtualDisplay* display)
                       void SetLuminosity(Int_t lumi = 1)
                       void SetBfield(Bool_t field = 1)
                       void SetSmearing(Bool_t val = 1)
                       void SetSmearMuonOpt(Int_t val = 3)
                       void SetSmearMuonFun(Int_t val = 2)
                       void SetSUSYcodeLSP(Int_t val = 66)
                       void SetTrackFinding(Bool_t val = 0)
                       void SetRun(Int_t run = 1)
                       void SetEvent(Int_t event = 1)
                       void SetMode(Int_t mode = 0)
                    TClass* Class()
                    TClass* IsA()
                       void ShowMembers(TMemberInspector& insp, char* parent)
                    ATLFast ATLFast(ATLFast&)


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.