5. Quick analysis ==================== 5.1 Analysis using qana.C ------------------------- If you would like to perform a quick analysis on your simulated data, there is also a handy tool available which simplifies a 'standard' analysis to a single line of code. The basis for this is the combiner class **PndSimpleCombiner**, the corresponding task **PndSimpleCombinerTask**, and the steering macro **qana.C**. Doing the analysis described above basically reduces to the command:: > root -l -b -q 'qana.C("signal_pid.root", 6.232, "J/psi->mu+ mu-;pbarpSystem->;J/psi pi+ pi-", "fit4c:fitvtx:mwin(J/psi)=0.8")' The result is the output root file **signal_pid_ana.root** with the contents:: root [1] .ls TFile** signal_pid_ana.root TFile* signal_pid_ana.root KEY: TFolder pnd;1 Main Folder KEY: TList BranchList;1 Doubly linked list KEY: TList TimeBasedBranchList;1 Doubly linked list KEY: FairFileHeader FileHeader;1 KEY: TTree ntp0;1 J/psi->mu+ mu- KEY: TTree ntp1;1 pbarpSystem->J/psi pi+ pi- KEY: TTree pndsim;1 /pnd_0 where ``ntp0`` and ``ntp1`` are the ntuples containing the reconstructed *J/psi* and *pbarpSystem* candidates. When running the macro without parameters (or with just "H" as first parameter) a (long) help text about the usage is printed:: qana.C -- ROOT macro for quick and simple analysis of PANDA simulation or data. USAGE: qana.C( , , , [parms], [nevt], [runnum], [mode], [runST] ) : input file name with PndPidCandidates; "help"/"H" shows analysis parameter options : pbar momentum; negative values are interpreted as -E_cm : decay pattern to be reconstructed, e.g. "phi -> K+ K-; D_s+ -> phi pi+"; decay file: auto generate from decay file [parms] : analysis options, e.g. "mwin=0.4:emin=0.1:pmin=0.1:qamc"; "auto": auto generate options; "qapart" runs particle QA [nevt] : number of events to analyze (default = 0 -> all) [runnum] : arbitrary integer run number (default = 0) [mode] : arbitrary integer mode number (default = 0) [runST] : run software trigger (default = false) Example 1 - Analyze EvtGen events : root -l -b -q 'qana.C("signal_pid.root", 6.23, "J/psi -> mu+ mu-; pbp -> J/psi pi+ pi-", "fit4c:mwin=0.6")' Example 2 - Analyze BOX generator events : root -l -b -q 'qana.C("box_pid.root", 15.0, "", "qapart")' available options for parameter [parms]; use colon separated setting 'par1=value1:par2=value2' ... - mwin : mass window for all composites - mwin(X) : mass window for composite 'X'; X (or cc of X) has to appear in decay - emin : minimum energy threshold for neutrals - pmin : minimum momentum threshold for charged - pid : PID criterion (al=All, vl=VeryLoose, lo=Loose, ti=Tight, vt=VeryTight, be=Best) for all species, e. g. pid=lo - pidX : PID criterion for specific species X = (e, mu, pi, k, p), e.g. pide=vl, pidmu=ti,.... - algo : PID algorithm (e.g. "algo=ideal,emc;drc;disc" or "algo=comb") for all species - algoX : PID algorithm for specific species X = (e, mu, pi, k, p) - fit4c[ K+ K-; D_s+ -> phi pi+", '!ntp0' would skip dump of TTree for phi->KK Available PID algorithms (FullSim / FastSim): - ideal = ideal particle id (PidAlgoIdealCharged / IdealPidProbability) - comb = combination of all following algos - mvd = Micro Vertex Detector (PidAlgoMvd / MvdPidProbability) - mdt = Muon Detector (PidAlgoMdtHardCuts / ScMdtPidBarrelProbability;ScMdtPidForwardProbability) - drc = DIRC detector (PidAlgoDrc / DrcBarrelProbability) - disc = Disc DIRC (PidAlgoDisc / DrcDiscProbability) - stt = Straw Tube Tracker (PidAlgoStt / SttPidProbability) - emc = EM-Calorimeter (PidAlgoEmcBayes / ScEmcPidFSProbability;ScEmcPidFwCapProbability;ScEmcPidBarrelProbability;ScEmcPidBwCapProbability) - rich = Forward RICH detector (PidAlgoRich / RichProbability) - scit = SciTil detector (PidAlgoSciT / n.a.) - ftof = Forward Time Of Flight (PidAlgoFtof / n.a.) The two most important parameters are the decay tree definition ``decay``, and the steering parameter string ``parms``, which are discussed a bit more in detail below. 5.1.1 The decay configuration parameter ++++++++++++++++++++++++++++++++++++++++++ Although the syntax of this parameter is quite intuitive, there are some rules which have to be followed for proper function. * The particle names used are those from *EvtGen*. * The decay chain has to be build up from generic particles: *e+, e-, mu+, mu+, pi+, pi-, K+, K-, p+, anti-p-, gamma*. * The sub-decays are separated by a semicolon. * The decay chain must be build bottom up, i.e. a particle used in a decay definition has to be defined before * *correct*: ``"pi0 -> gamma gamma; omega -> pi+ pi- pi0"`` * *wrong*: ``"omega -> pi+ pi- pi0; pi0 -> gamma gamma"`` * As default, charge conjugation is applied automatically. If this is not wanted, the key word **nocc** can be added at the end of each decay. This can be useful, if particle and anti-particle go to different decay channels. * ``"D0 -> K- pi+; pbarpSystem -> D0 anti-D0"`` reconstructs both *D0 -> K- pi+* and *anti-D0 -> K+ pi-* and combines them to a *pbar-p-system*. Only one n-tuple is stored for ``D0`` and ``anti-D0``, distinguishable by the branch ``xpdg`` * ``"D0 -> K- pi- nocc; pbarpSystem -> D0 anti-D0"`` fails to run, since the list ``anti-D0`` is missing. * ``"D0 -> K- pi+ nocc; anti-D0 -> pi+ pi- nocc; pbarpSystem -> D0 anti-D0"`` reconstruct ``pbarpSystem`` from ``D0`` and ``anti-D0`` from different decay channels. * For each sub-decay an n-tuple is stored by default. These are named ntp0, ntp1, ntp2,... with the index according to the position of the sub-decay in the decay configuration string. This behaviour can be modified in the ``parms`` setting. 5.1.2 The steering parameter setting ++++++++++++++++++++++++++++++++++++++ There is a default behavior how the analysis is carried out, which basically is restricted to pure combinatorics. Many aspects like fitting, PID, mass windows, etc, can be adjusted with the *parms* parameter. Some of the settings are interpreted by the macro *quickana.C*, some by the combiner class *PndSimpleCombiner*, and the others by the task class *PndSimpleCombinerTask*. The following table lists all available options. **The individual parameters are separated by a colon (:).** .. list-table:: Steering parameters :header-rows:1 * - Parameter - Meaning * - fit4c[x GeV for photon candidates. * - pmin=x - Minimum momentum cut p>x GeV/c for tracks. * - pid=x - Common PID criterion for all particle species. x ∈ [All, VeryLoose, Loose, Tight, VeryTight, Best]. Default is All. * - pidY=x - PID criterion for species Y. x ∈ [All, VeryLoose, Loose, Tight, VeryTight, Best], Y ∈ [e, mu, pi, k, p]. Examples: pidmu=VeryTight, pidk=Loose. * - algo=x - Common PID algorithm configuration x, e.g algo PidAlgoEmcBayes;PidAlgoDrc. Default for all species is PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts;PidAlgoSciT;PidAlgoRich * - algoY=x - PID algorithm setting for species Y, Y ∈ [e, mu, pi, k, p]. Example\: algoe=PidAlgoEmcBayes:algok=PidAlgoDrc;PidAlgoDisc;PidAlgoStt. * - gtrk=x - GoodTrackSelector mode. x ∈ [All, VeryLoose, Loose, Tight, VeryTight] * - gphot=x - GoodPhotonSelector mode. x ∈ [All, VeryLoose, Loose, Tight, VeryTight] * - ebrem - Use the electron lists with bremsstrahlung correction (ElectronBrem). * - qamc - Stores the MC truth list in the n-tuple nmc. * - qaevtshape - Stores event shape variables in *each* n-tuple. * - nevt - Stores *additional n-tuple with event (shape) variables*; one entry per event. * - qapart - Runs the task **PndParticleQATask**, useful for single particle studies. Creates n-tuples *ntp* (charged), *ntpn* (neutrals), *nmc* (MC truth list) * - !ntpx - Prevents storing the n-tuple from the x-th decay definition. Example: !ntp0:!ntp1 suppresses storage of the two first decay definitions in decay chain. * - !mc - Prevents **PndParticleQATask** to store the MC information in n-tuple nmc. * - !neut - Prevents **PndParticleQATask** to store the single neutral information in n-tuple ntpn. * - !chrg - Prevents **PndParticleQATask** to store the charged track information in n-tuple ntp. 5.2 On-the-fly Fast Simulation and Analysis using qfsana.C ------------------------------------------------------------ In case your want to take a really quick look on how a particular decay channel reconstructs in terms of acceptance or kinematic distributions, you can use make use of Fast Simulations together with the above mentioned simplified combinatorics by using the macro **qfsana.C**. Instead of simulated data, you just need a EvtGen decay file, or, for background studies not even that. The complete simulation and analysis described above for signal and background can be performed with two lines of code:: > root -l -b -q 'qfsana.C("signal", // prefix "pp_jpsi2pi_jpsi_mumu.dec", // decay file 6.232, // pbar momentum "J/psi->mu+ mu-;pbarpSystem->J/psi pi+ pi-", // decay pattern 1000, // number of events "fit4c:fitvtx:mwin=0.8")' // config parameters > root -l -b -q 'qfsana.C("bkg", "DPM", // other generator 6.232, "J/psi->mu+ mu-;pbarpSystem->J/psi pi+ pi-", 10000, "fit4c:fitvtx:mwin=0.8")' The resultant files contains the same n-tuples as shown above. Documentation about the names of the branches can be found under the description of PndRhoTupleQA (https://panda-wiki.gsi.de/foswiki/bin/view/Computing/PandaRootAnalysisJuly13#PndRhoTupleQA). The parameters for this macro are slightly different than those for **qana.C**, however, the parameters for *decay* and *parms* are the same as above:: qfsana.C -- ROOT macro for quick and simple fast simulation and analysis for PANDA USAGE: qfsana.C( , , [mom], [nevt], [decay], [parms], [runST], [runnum], [mode] ) : output file names prefix; "help"/"H" shows analysis parameter options : EvtGen decfile "xxx.dec" or "xxx.dec:iniResonance"; DPM/FTF/BOX uses DPM, FTF or BOX generator - DPM settings: DPM = inelastic only, DPM1 = inel. + elastic, DPM2 = elastic only - FTF settings: FTF = inelastic only, FTF1 = inel. + elastic - BOX settings: "BOX:type[pdgcode,mult]"; optional range/value with "p/tht/ctht/phi[min,max]/[value]" separated by ":", e.g. "BOX:type[211,2]:p[1,5]:tht[45]:phi[90,210]" [mom] : pbar momentum (negative values are interpreted as -E_cm) for EvtGen/DPM/FTF; BOX generator: maximum particle momentum (if not specified above) (default = 15.0) [nevt] : number of events to simulate and analyze (default = 1000) [decay] : decay pattern to be reconstructed, e.g. "phi -> K+ K-; D_s+ -> phi pi+"; "auto"/decfile: generate from decay file; "": fastsim only w/o reco (default = auto) [parms] : analysis options, e.g. "mwin=0.4:emin=0.1:pmin=0.1:qamc"; "auto": auto generate options; "qapart" runs particle QA; "persist" saves PndPidCandidates (default = auto) [runnum] : arbitrary integer run number (default = 0) [mode] : arbitrary integer mode number (default = 0) [runST] : run software trigger (default = false) Example 1 - EvtGen sim/analysis : root -l -b -q 'qfsana.C("jpsi2pi", "pp_jpsi2pi_jpsi_mumu.dec:pbarpSystem", 6.23, 1000, "J/psi -> mu+ mu-; pbp -> J/psi pi+ pi-", "fit4c:mwin=0.6")' Example 2 - EvtGen auto mode : root -l -b -q 'qfsana.C("jpsi2pi_auto", "pp_jpsi2pi_jpsi_mumu.dec", 6.23)' Example 3 - Particle QA BOX gen : root -l -b -q 'qfsana.C("single_kplus", "BOX:type[321,1]", 8.0, 1000, "", "qapart")' Example 4 - Simulation output DPM : root -l -b -q 'qfsana.C("bkg", "DPM", 6.23, 1000, "", "")'