.. _ht: .. currentmodule:: pyrecola Higgs Triplet Extension of the Standard Model ============================================= This Higgs Triplet Extension of the Standard Model (*HT*) is a simple extended Higgs sector with one additional real vevless SU(2) Higgs triplet field :m:`\Sigma` with hypercharge :math:`Y_\Sigma=0`. Our conventions follow mainly :cite:`FileviezPerez:2008bj`. See also :cite:`Chen:2008jg`. In the vevless case the most general renormalizable scalar potential reads .. math:: V_{\mathrm{HT}} = -m_\Phi^2 \Phi^\dagger \Phi -m_\Sigma^2 F +\frac{\lambda}{4} \left(\Phi^\dagger \Phi\right)^2 +\frac{b_4}{4} F^2 +\frac{a_2}{2} \Phi^\dagger \Phi F, :label: eq:htpot with :math:`\Phi` being the *SM* Higgs doublet and :math:`F` being defined as .. math:: F = 2 \mathrm{Tr} \left[ \Sigma\cdot\Sigma\right]. Parameter conventions ********************* Parameters of the model can be set with the generic function :py:meth:`set_parameter_rcl`. All parameters in the HT are defined to be real and we choose the following set of physical parameters: .. list-table:: :widths: 6 10 4 :header-rows: 1 :align: center * - Basis - HT potential - Gauge part * - before SSB - :m:`m_\Phi`, :m:`m_\Sigma`, :m:`\lambda`, :m:`a_2`, :m:`b_4` - :m:`g`, :m:`g^\prime` * - Recola2 input - :m:`M_{\mathrm{H}_1}`, :m:`M_{\mathrm{H}^\pm}`, :m:`a_2`, :m:`b_4`, :m:`M_\mathrm{W}` - :m:`\alpha_\mathrm{em}`, :m:`M_\mathrm{Z}` The default values are .. list-table:: :widths: 6 12 8 :header-rows: 1 :align: center * - Parameter - |recola| identifier - default value * - :m:`M_{\mathrm{H}_1}` - ``'MH1'`` - 125 * - :m:`M_{\mathrm{H}^\pm}` - ``'MHP'`` - 500 * - :m:`a_{2}` - ``'a2'`` - 0.1 * - :m:`b_{4}` - ``'b4'`` - 0.2 The fields extend the ones in the :ref:`SM ` by .. list-table:: :widths: 10 12 :header-rows: 1 :align: center * - Fields - |recola| identifier * - :m:`H_\mathrm{1}` - :code:`'H1'` * - :m:`H_\mathrm{2}` - :code:`'H2'` * - :m:`H^+` - :code:`'H+'` * - :m:`H^-` - :code:`'H-'` where :m:`H_\mathrm{1}` is the *lighter* Higgs-boson which typically takes the role of the SM one. Power counting and renormalization ********************************** The model has been implemented with a power counting (see :ref:`SM power counting `) that assumes the couplings :m:`a_2`, :m:`b_2`` to scale as :m:`\mathrm{QED}^2`. The renormalization has been performed in the complete on-shell scheme with the possibility to switch between standard EW renormalization schemes. The new Higgs-potential parameters are renormalized MSbar. Concerning the renormalization see also :cite:`Blank:1997qa`. The mass correction to the second neutral Higgs-boson can be printed by increasing the print level :py:meth:`set_print_level_parameters_rcl`. Snippet code using the HT ************************* .. tabs:: .. code-tab:: py from pyrecola import * set_output_file_rcl('*') set_print_level_squared_amplitude_rcl(2) # Change HT parameters set_parameter_rcl("MH1", 125.) set_parameter_rcl("MHP", 350.) set_parameter_rcl("a2", 0.3) set_parameter_rcl("b4", -0.7) # enable to draw off-shell currents # set_draw_level_branches_rcl(1) define_process_rcl(1, 'e+ e- -> H+ H-', 'NLO') generate_processes_rcl() p1 = [500., 0., 0., 500.] p2 = [500., 0., 0., -500.] # generate a sample PSP using RAMBO p = set_outgoing_momenta_rcl(1, [p1, p2]) # compute tree squared and tree one-loop interference compute_process_rcl(1, p, 'NLO') # get all different contributions (pow=[n,m,o] == gs^n e^m k^o) A0 = get_squared_amplitude_rcl(1, 'LO', pow=[0, 4]) # born-virtual interference (EW corrections) A1 = get_squared_amplitude_rcl(1, 'NLO', pow=[0, 6]) print("A0:", A0) print("A1:", A1) reset_recola_rcl() .. code-tab:: fortran program main use recola implicit none integer, parameter :: dp = kind (23d0) real(dp) :: p(0:3,1:4), A0, A1 call set_parameter_rcl("MH1", complex(125d0,0d0)) call set_parameter_rcl("MHP", complex(350,0d0)) call set_parameter_rcl("a2", complex(0.3d0,0d0)) call set_parameter_rcl("b4", complex(-0.7d0,0d0)) call set_output_file_rcl('*') call set_print_level_squared_amplitude_rcl(2) ! enable to draw off-shell currents ! call set_draw_level_branches_rcl(1) call define_process_rcl(1, 'e+ e- -> H+ H-', 'NLO') call generate_processes_rcl p(:,1) = [500d0, 0d0, 0d0, 500d0] p(:,2) = [500d0, 0d0, 0d0, -500d0] ! generate a sample PSP using RAMBO call set_outgoing_momenta_rcl(1, p(:,1:2), p) ! compute tree squared and tree one-loop interference call compute_process_rcl(1, p, 'NLO') call get_squared_amplitude_rcl(1,[0,4], 'LO' , A0) call get_squared_amplitude_rcl(1,[0,6], 'NLO', A1) write(*,*) "A0:", A0 write(*,*) "A1:", A1 end program main .. code-tab:: c++ #include "recola.hpp" #include int main(int argc, char *argv[]) { Recola::set_output_file_rcl("*"); Recola::set_print_level_squared_amplitude_rcl(2); Recola::set_parameter_rcl("MH1", 125.); Recola::set_parameter_rcl("MHP", 350.); Recola::set_parameter_rcl("a2", 0.3); Recola::set_parameter_rcl("b4", -0.7); // enable to draw off-shell currents // Recola::set_draw_level_branches_rcl(1); Recola::define_process_rcl(1, "e+ e- -> H+ H-", "NLO"); // generate it Recola::generate_processes_rcl(); // generate a sample PSP using RAMBO double pin[2][4] = {{500., 0., 0., 500.}, {500., 0., 0., -500.}}; double p[4][4]; Recola::set_outgoing_momenta_rcl(1, pin, p); // compute tree squared and tree one-loop interference double A2[2]; Recola::compute_process_rcl(1, p, "NLO", A2); double A0,A1; int pow[2] = {0, 4}; Recola::get_squared_amplitude_rcl(1, pow, "LO", A0); pow[1] = 6; Recola::get_squared_amplitude_rcl(1, pow, "NLO", A1); std::cout << "A0: " << A0 << std::endl; std::cout << "A1: " << A1 << std::endl; return 0; } Releases ^^^^^^^^ * :ht:`2.2.3` **newest** * :ht:`2.2.2` UFO model files ^^^^^^^^^^^^^^^ * :ht:`UFO` .. rubric :: References .. bibliography:: ../references.bib :filter: docname in docnames