Changes between Version 3 and Version 4 of LOEventGenerationBias


Ignore:
Timestamp:
Aug 20, 2016, 4:02:19 AM (8 years ago)
Author:
Valentin Hirschi
Comment:

--

Legend:

Unmodified
Added
Removed
Modified
  • LOEventGenerationBias

    v3 v4  
    88The motivation for using a bias function typically falls in one of the following two categories:
    99 * a) Producing smoother distributions for the tail of a particular observable. This means that physical results obtained in presence of the bias will be identical but sampled differently.
    10  * b) One wants to modify the integrand so as to really impact the physical results. This can be useful for a plethora of applications: ad-hoc unitarisation of the matrix elements, merging weights, inclusion of higher order contributions, etc..
     10 * b) One wants to modify the integrand so as to really impact the physical results. This can be useful for a plethora of applications: ad-hoc unitarisation of the matrix elements, merging weights, inclusion of higher order contributions, implementation of customized cuts, on-the-flight plotting, etc..
    1111
    1212Before going into the detail of the usage of the bias module and the instructions for building your customized bias, we present here an example of their use for the case a) above to smoothen the distribution of the jet transverse momentum for the process 'p p > e+ e- j'.
     
    2727
    2828== Usage of a bias module
     29
     30In order to activate a given bias module for a MadEvent run, you must specify it in the run card in the entry '{{{bias_module}}}'.
     31Alternatively one can use the syntax '{{{set bias_module <chosen_module>}}}' to write this entry when MG5aMC offers you to edit the cards.
     32
     33The chosen bias module can either be one of the bias modules natively shipped with MG5aMC ( like '{{{ptj_bias}}}' ) or an absolute path pointing to a folder that conforms to the format required to implement your own custom bias function (see next section). In any case, in order to turn off any bias, you can set the '{{{bias_module}}}' run_card entry to '{{{None}}}'.
     34
     35Similarly, the bias parameters (all of double precision type) can be specified with a dictionary-like syntax in the run_card entry '{{{bias_parameters}}}'. For example:
     36
     37{{{ {'ptj_bias_target_ptj': 300.0, 'ptj_bias_enhancement_power': 4.0} = bias_parameters }}}
     38
     39Alternatively, the commands '{{{set bias_parameter <param_dictionary>}}}' and '{{{add bias_parameter <param_name> <value>}}}' can be used when MG5aMC offers you to edit the cards.
     40
     41Finally, notice that when you choose to use a bias module that does *not* impact the physical results, the unweighted events produced will be re-weighted so as to undo the effect of the bias.
     42As a result the events produced by MadEvent will weighted in this case (obviously unweighting them again would return the original unweighted distribution would would have obtained without the use of a bias in the first place). 
     43
     44== Writing your own bias module
     45
     46In order to develop your custom bias module, you must create a directory that bears the name of the module. For the sake of simplicity, we will illustrate the requirements for a custom bias module with the example of a module that applies on each event a CKKW-L type of merging weight obtained from Pythia8. We will call this bias module '{{{PY8_CKKW}}}'.
     47You can also have a look at the structure of the much simpler implementation of the '{{{ptj_bias}}}' module in '{{{<MG5aMC_distrib_path>/Template/LO/Source/BIAS/ptj_bias}}}'.
     48
     49The directory '{{{PY8_CKKW}}}' must then minimally contain the two files '{{{PY8_CKKW/PY8_CKKW.f}}}' and '{{{PY8_CKKW/makefile}}}'. The optional file '{{{PY8_CKKW/bias_dependencies}}}' can also be present.
     50On top of these files, the user is of course free to include any additional resource that its module will use, but these will not play a role in the interface of the module to MG5_aMC.
     51
     52=== Mandatory file : '{{{<BIAS_MODULE_NAME>/<BIAS_MODULE_NAME>.f}}}'
     53
     54As you can see, this fortran file must bear the name of the module with the suffix '{{{.f}}}' and contain the subroutine '{{{bias_wgt}}}' with the prototype
     55
     56{{{
     57subroutine bias_wgt(p, original_weight, bias_weight)
     58include '../../nexternal.inc'
     59          double precision p(0:3,nexternal)
     60          double precision original_weight, bias_weight
     61}}}
     62
     63Where the array '{{{p}}}' specifies the kinematics of the events, and the '{{{original_weight}}}' is also provided, in case your bias function was not to be multiplicative.
     64The last argument '{{{bias_weight}}}' is intended as the output.
     65
     66The parameters of your bias module can be specified by including the following comment in the body of the subroutine:
     67
     68{{{
     69C  parameters = {'OneParam': 0.0,
     70C                           'AnotherParam': 0.0,
     71C                           'YetAnotherOne':2.0}
     72}}}
     73
     74You must makes sure that your subroutine defines double-precision variables of the same name in its body.
     75The value of these parameters chosen by the user in the run_card can be accessed by including the following include file (after all variable declarations):
     76
     77{{{
     78          include '../bias.inc'
     79}}}
     80
     81Finally the module must specify the following common block:
     82
     83{{{
     84          double precision stored_bias_weight
     85          data stored_bias_weight/1.0d0/
     86          logical impact_xsec, requires_full_event_info
     87C        Decides whether the bias must be corrected for after unweighting
     88          data impact_xsec/.True./
     89C        Does your bias need the full information about the event
     90C         (color, resonances, helicities, etc..)
     91          data requires_full_event_info/.True./
     92          common/bias/stored_bias_weight,impact_xsec,
     93     &                requires_full_event_info
     94}}}
     95
     96The initialization of the logicals '{{{impact_xsec}}}' and '{{{requires_full_event_info}}}' are very important:
     97
     98 * If '{{{impact_xsec}}}' is set to '{{{.False.}}}', MG5aMC will make sure to reweight the unweighted events produced in presence of the bias so as to undo the physical impact of the bias function (without affecting the resulting event distribution of course). On the contrary, if you want the effect of your bias to persist on the cross-section and differential distribution, you must set this parameter to '{{{.True.}}}'.
     99 * If '{{{requires_full_event_info}}}' is set to '{{{.True.}}}' will make sure to fully process all events sent to your bias module so that color, helicity, resonances information are available. If your module only relies on the kinematics of the event, you should set this parameter to '{{{.False.}}}' so as to speed up the code.
     100
     101Finally, you should store your desired bias once computed in the subroutine argument '{{{bias_weight}}}'.
     102
     103Your module can have access to much more information that just the kinematics specified in the subroutine argument. To access it, simply add the following include file
     104
     105{{{
     106'          include '../../lhe_event_infos.inc'
     107}}}
     108
     109There are also other include files that you can add to access other information like the cuts values, the run card parameters, etc...
     110
     111The illustrative example file {{{PY8_CKKW/PY8_CKKW.f}}} highlights the use of this additional information. For completion we also show the file '{{{py8_mg5amc_bias_interface.cc}}}' which implements the C function '{{{py8_bias_weight}}}' used in '{{{PY8_CKKW.f}}}'.
     112
     113=== Mandatory file : '{{{<BIAS_MODULE_NAME>/makefile}}}'
     114
     115This makefile will be called by MG5aMC to compile your module into a library when the user choses to use it.
     116It must comply with the following rules:
     117
     118 * a) Include the general makefile options of MG5aMC with:
     119
     120{{{
     121BIASDIR   = ../
     122SOURCEDIR = ../..
     123include $(SOURCEDIR)/make_opts
     124}}}
     125
     126 * b) Define the following target '{{{clean:}}}':
     127
     128{{{
     129clean:
     130        $(RM) *.o $(BIASLIBDIR)$(BIASLIBRARY)
     131}}}
     132
     133 * c) Define the target '{{{<BIAS_MODULE_NAME>:}}}' that creates the library '{{{libbias.a}}}' as follows:
     134
     135{{{
     136<BIAS_MODULE_NAME>: <your list of object.o files>
     137        $(call CREATELIB, $(BIASLIBDIR)$(BIASLIBRARY), $^)
     138}}}
     139
     140  * d) Define the target '{{{requirements:}}}' that prints out 'VALID' if the requirements of your module are met. For example
     141
     142{{{
     143requirements:   
     144ifeq ($(shell $(call CHECK_MG5AMC_VERSION,2.4.2)),True)
     145ifeq ($(PYTHIA8_PATH),NotInstalled)
     146        @echo "Error: This module requires Pythia8. Consider installing it with the command 'install pythia8'."
     147else
     148        @echo "VALID"
     149endif
     150else
     151        @echo "Error:: MG5aMC is not recent enough (found "$(MG5AMC_VERSION)")"
     152        @echo "FAIL"
     153endif
     154}}}
     155
     156which uses the macro '{{{CHECK_MG5AMC_VERSION}}}' defined in '{{{../make_opts}}}'.
     157
     158You can then add to this makefile whatever target and dependencies that are necessary for compiling your module.
     159
     160We show here the full makefile corresponding to our example module '{{{PY8_CKKW}}}'.
     161
     162=== Optional file : '{{{<BIAS_MODULE_NAME>/bias_dependencies}}}'
     163
     164Finally your bias module might link against additional libraries that the MadEvent fortran executables should also link to when using the module.
     165To specify them, you can add the file '{{{bias_dependencies}}}' to your bias module directory, where you define the variable 'BIASDEPENDENCIES' which will be automatically added to MadEvent makefiles.
     166
     167Here is what this optional file '{{{bias_dependencies}}}' looks like in our example module '{{{PY8_CKKW}}}':
     168
     169{{{
     170PY8LIBS     = $(shell $(PYTHIA8_PATH)/bin/pythia8-config --libs)
     171PY8INCLUDES = $(shell $(PYTHIA8_PATH)/bin/pythia8-config --includedir)
     172PY8CXXFLAGS = $(shell $(PYTHIA8_PATH)/bin/pythia8-config --cxxflags)
     173
     174BIASDEPENDENCIES = $(PY8LIBS) $(PY8CXXFLAGS)
     175}}}
     176
     177Note that the environment variable '{{{PYTHIA8_PATH}}}' will be in this case provided by '{{{make_opts}}}', also sourced in MadEvent makefiles.