Skip to content
Snippets Groups Projects

Starterkit: Update line writing tutorial; Add prototype Bs -> Jpsi phi lines

Merged Marian Stahl requested to merge mstahl_starterkit into master
Files
10
###############################################################################
# (c) Copyright 2019 CERN for the benefit of the LHCb Collaboration #
# #
# This software is distributed under the terms of the GNU General Public #
# Licence version 3 (GPL Version 3), copied verbatim in the file "COPYING". #
# #
# In applying this licence, CERN does not waive the privileges and immunities #
# granted to it by virtue of its status as an Intergovernmental Organization #
# or submit itself to any jurisdiction. #
###############################################################################
"""Define HLT2 line for ``Lambda_b0 -> Lambda_c+ pi+``.
With ``Lambda_c+ -> p+ K- pi+``.
An example input file and Moore configuration is also given at the bottom of
this file, so that it can be run as-is.
"""
import Functors as F
from Functors.math import in_range
from GaudiKernel.SystemOfUnits import GeV, MeV, mm
from Moore.config import register_line_builder
from Moore.lines import Hlt2Line
from RecoConf.reconstruction_objects import (
make_pvs_v2 as make_pvs,
upfront_reconstruction,
)
# For code inside files under Hlt2Conf/python/lines you should reference these
# two modules using relative imports:
#
# from ..standard_particles import make_has_rich_long_kaons
from Hlt2Conf.standard_particles import (
make_has_rich_long_kaons,
make_has_rich_long_pions,
make_has_rich_long_protons,
)
from Hlt2Conf.algorithms_thor import ParticleCombiner, ParticleFilter, require_all
all_lines = {}
def filter_protons(particles, pvs, pt_min=0.5 * GeV, mipchi2_min=9,
dllp_min=5):
cut = require_all(
F.PT > pt_min,
F.MINIPCHI2(pvs) > mipchi2_min,
F.PID_P > dllp_min,
)
return ParticleFilter(particles, F.FILTER(cut))
def filter_kaons(particles, pvs, pt_min=0.5 * GeV, mipchi2_min=9, dllk_min=5):
cut = require_all(
F.PT > pt_min,
F.MINIPCHI2(pvs) > mipchi2_min,
F.PID_K > dllk_min,
)
return ParticleFilter(particles, F.FILTER(cut))
def filter_pions(particles, pvs, pt_min=0.5 * GeV, mipchi2_min=9, dllk_max=5):
cut = require_all(
F.PT > pt_min,
F.MINIPCHI2(pvs) > mipchi2_min,
F.PID_K < dllk_max,
)
return ParticleFilter(particles, F.FILTER(cut))
def make_lambdacs(protons,
kaons,
pions,
pvs,
two_body_comb_maxdocachi2=9.0,
comb_m_min=2080 * MeV,
comb_m_max=2480 * MeV,
comb_pt_min=2000 * MeV,
comb_maxdoca=0.1 * mm,
vchi2pdof_max=10,
bpvvdchi2_min=25):
two_body_combination_code = F.MAXDOCACHI2CUT(two_body_comb_maxdocachi2)
combination_code = require_all(
in_range(comb_m_min, F.MASS, comb_m_max),
F.SUM(F.PT) > comb_pt_min,
F.MAXDOCACUT(comb_maxdoca),
)
vertex_code = require_all(
F.CHI2DOF < vchi2pdof_max,
F.BPVFDCHI2(pvs) > bpvvdchi2_min,
)
return ParticleCombiner(
[protons, kaons, pions],
DecayDescriptor="[Lambda_c+ -> p+ K- pi+]cc",
Combination12Cut=two_body_combination_code,
CombinationCut=combination_code,
CompositeCut=vertex_code,
)
def make_lambdabs(lcs,
pions,
pvs,
comb_m_min=5000 * MeV,
comb_m_max=7000 * MeV,
comb_pt_min=4000 * MeV,
comb_maxdoca=0.1 * mm,
vchi2pdof_max=10,
bpvvdchi2_min=25):
combination_code = require_all(
in_range(comb_m_min, F.MASS, comb_m_max),
F.SUM(F.PT) > comb_pt_min,
F.MAXDOCACUT(comb_maxdoca),
)
vertex_code = require_all(
F.CHI2DOF < vchi2pdof_max,
F.BPVFDCHI2(pvs) > bpvvdchi2_min,
)
return ParticleCombiner(
[lcs, pions],
DecayDescriptor="[Lambda_b0 -> Lambda_c+ pi-]cc",
CombinationCut=combination_code,
CompositeCut=vertex_code,
)
@register_line_builder(all_lines)
def lbtolcpi_lctopkpi_line(name="Hlt2LbToLcpPim_LcToPpKmPipLine", prescale=1):
pvs = make_pvs()
protons = filter_protons(make_has_rich_long_protons(), pvs)
kaons = filter_kaons(make_has_rich_long_kaons(), pvs)
pions = filter_pions(make_has_rich_long_pions(), pvs)
lcs = make_lambdacs(protons, kaons, pions, pvs)
lbs = make_lambdabs(lcs, pions, pvs)
return Hlt2Line(
name=name,
algs=upfront_reconstruction() + [lbs],
prescale=prescale,
)
# Moore configuration
from Moore import options, run_moore
# In a normal options file, we would import the line from Hlt2Conf where it is
# defined
# from Hlt2Conf.lines.LbToLcPi import lbtolcpi_lctopkpi_line
# Temporary workaround for TrackStateProvider
from RecoConf.global_tools import stateProvider_with_simplified_geom
public_tools = [stateProvider_with_simplified_geom()]
def all_lines():
return [lbtolcpi_lctopkpi_line()]
options.set_input_and_conds_from_testfiledb('Upgrade_MinBias_LDST')
options.input_raw_format = 4.3
options.evt_max = 100
options.control_flow_file = 'control_flow.gv'
options.data_flow_file = 'data_flow.gv'
run_moore(options, all_lines, public_tools)
Loading