Multithreading in generators
Goal of this page
CMS has implemented a variety of multithreading methods in the event generation step under the CMSSW workflow. After reading this page, you will
- Learn how to configure the fragment or the python configuration file to enable these settings.
- Understand how these features actually function to help speed up the program and increase CPU efficiency.
Contents
Outline of event generation in CMS
Basic review
This section provides a review of event generation in CMS. Some aspects will be useful to understand the multithreading implementation described later. Details of event generation are documented in
WorkBookGeneration.
The generators can be classified as the
general-purpose generators (e.g. Pythia8, Herwig7, ...) or
matrix element (ME) generators (MG5_aMCatNLO, Powheg, ...). There are two ways in CMS to generate GEN events. In the first method, a single general-purpose generator is used to directly produce GEN events. On the second case, the ME step and the event hadronization are separated. The ME step aims to produce LHE-format events, which either comes from an existing LHE file, or from the ME calculation by a specific ME generator. The LHE events are then interfaced to the general-purpose generator for hadronization.
In the CMS workflow, the two types of generators are treated differently. The ME generators always run in a step external to CMSSW, thus the LHE events production is not controlled by the
cmsRun
process. The general-purpose ones, as a contrary, has a dedicated interface to the CMSSW framework, therefore the generation of each event is scheduled by CMSSW.
Below is a review on how to set the configurations for the two cases.
In the first case when we only need one general-purpose generator to produce events, GEN configuration is specified in one module. A frequently-used example (at the time of writing) is
Pythia8GeneratorFilter
that uses Pythia8 generator to directly produces GEN events. Other examples may include e.g.
Pythia8PtGun
that make use of specific particle gun. All such modules are based on a
EDFilter
class, thus event-generation is under the control of
cmsRun
. Commands passing to the generator to produce a given physics process are declared in the module. Here is one configuration example:
QCD_Pt-15To7000_TuneCUETP8M1_Flat_14TeV-pythia8_cff.py
.
In the second case, if the LHE events are not from existing LHE records but from a dedicated generator, we need the ME generators and the general-propose generator to cooperate. The first step, i.e. the LHE event generation, is always triggered by the
ExternalLHEProducer
module. The module can be found in the fragment or in the configuration file, always appearing as the below format.
process.externalLHEProducer = cms.EDProducer("ExternalLHEProducer",
args = cms.vstring('/cvmfs/path/to/tarball.tar.xz'),
nEvents = cms.untracked.uint32(10000),
numberOfParameters = cms.uint32(1),
outputFile = cms.string('cmsgrid_final.lhe'),
scriptName = cms.FileInPath('GeneratorInterface/LHEInterface/data/run_generic_tarball_cvmfs.sh'),
)
The configuration of
ExternalLHEProducer
includes a gridpack. The gridpack is essentially a pre-configured, compiled and sealed generator, containing necessary phase-space information for a given process, which helps to produce LHE events out of the box. Practically, it is doable to unzip the gridpack manually and produce LHE events directly from it. The
ExternalLHEProducer
module will execute a shell script given by
scriptName
(which is always unchanged), latter will first unzip the gridpack locally, then start to produce LHE events out from it. After the script has completed and the LHE file is created, the module then transforms the LHE file into the EDM-format ROOT file. Here is one configuration example:
WTolNu01234Jets_5f_LO_MLM_Madgraph_LHE_13TeV_cff.py
.
The EDM-format LHE events are then interfaced to a general-propose generator (e.g. Pythia8 or Herwig7) to produce GEN events, by using modules such as
Pythia8HadronizerFilter
. Similar to
Pythia8GeneratorFilter
, it is also based on
EDFilter
as part of the CMSSW infrastructure. See here as a configuration example:
Hadronizer_TuneCP5_13TeV_MLM_5f_max4j_LHE_pythia8_cff.py
.
Configuration for threading
In this section, we will briefly introduce the thread configuration and how it works in CMSSW.
The configuration for threading is provided by the following line
process.options.numberOfThreads=cms.untracked.uint32(4)
process.options.numberOfStreams=cms.untracked.uint32(0)
where you ask 4 threads to process the events, and the number of streams (if given 0) will be set to the number of threads, enabling the multithreading workflow with 4 streams. You can also append
--Threads 4
to the
cmsDriver.py
utility, so the above block will be generated in the python configuration file.
The above setting specifies the number of streams to be activated. The implementation details on multithreading will, practically, be given to the module to decide. CMSSW implements several modes for multithreading, each of them separately derives the base class
EDProducer
,
EDFilter
, and
EDAnalyzer
for different types of event processing. For the threaded frameworks, you can check
MultithreadedFrameworkDesignDiscussions for more details, and we will come back to those concepts when necessary in what follows.
As the
ExternalLHEProducer
module processes events in a way external to CMSSW, it passes the number of threads N to the launch script the module calls in the runtime, which then passes the number of threads to the script
runcmsgrid.sh
inside the gridpack, leaving the gridpack to decide whether and how multithreading will be activated according to the number of threads.
The following techniques are essentially the enhancement to multithreading in a number of ways, designed for the GEN step. The goal is to help further accelerate the speed and improve the CPU efficiency.
Please make sure to work on CMSSW with version >= CMSSW_10_6_18.
Multi-thread methods
Concurrent ExternalLHEProducer
The concurrent ExternalLHEProducer technique provides a general multi-thread solution to the
ExternalLHEProducer
module.
For
ExternalLHEProducer
, the feature is activated by setting its parameter
generateConcurrently
as True (the default is False). The parameter is not set by default, so you may need to add the parameter manually, as the following shows:
process.externalLHEProducer = cms.EDProducer("ExternalLHEProducer",
...
scriptName = cms.FileInPath('GeneratorInterface/LHEInterface/data/run_generic_tarball_cvmfs.sh'),
generateConcurrently = cms.untracked.bool(True),
)
You can also apply the setting by adding the following argument to
cmsDriver.py
, while using the original fragment.
--customise_commands='process.externalLHEProducer.generateConcurrently = cms.untracked.bool(True)'
Underlying concept
By setting
generateConcurrently
to True, the code will execute the script in multiple separate instances at the same time, each as an individual process and does not communicate with each other. The multithreading is implemented by the TBB library. The number of threads N is as specified in the configuration file. The total events are divided by N and assigned to each thread. The number of threads is then reset as 1 in each instance, to ensure that each individual instance will process events in a single-thread mode.
Use case
The concurrent ExternalLHEProducer feature is applicable to all
ExternalLHEProducer
module to produce LHE events. For those using a more light-weighted gridpacks, the multiple processes launched in parallel will not cause a heavy burden on the I/O or RAM resources. For such case, the multithreading feature may significantly reduce the event-generation time.
For some use cases, the gridpack can be very large in size or requires large RAM when processing (for instance, the MG5_aMCatNLO gridpack with multijet final states in the process / with a complex NLO process). The concurrent feature may not have expected advantage in these cases.
Comparison of performances across all multithreading methods is always recommended before launching the production. The last section performs some case studies as a reference.
Internal multithreading for MG5_aMCatNLO at LO
If the
ExternalLHEProducer
module uses an MG5_aMCatNLO gridpack in the LO mode, you can simply set the internal multithreading mode developed by the MG5_aMCatNLO team. CMS provides an interface to this feature in a specific launch script. To activate the feature, replace the default
scriptName
by the script shown as follows.
process.externalLHEProducer = cms.EDProducer("ExternalLHEProducer",
...
scriptName = cms.FileInPath('GeneratorInterface/LHEInterface/data/run_generic_tarball_cvmfs_madgraphLO_multithread.sh'),
)
You can also apply the setting by adding the following argument to
cmsDriver.py
, while using the original fragment.
--customise_commands="process.externalLHEProducer.scriptName = cms.FileInPath('GeneratorInterface/LHEInterface/data/run_generic_tarball_cvmfs_madgraphLO_multithread.sh')"
MG5_aMCatNLO gridpack at NLO, on the other hand, do not require any special configuration to enable the multi-thread.
Underlying concept
By default, MG5_aMCatNLO supports multithreading by setting the
nb_core
. This is how MG5_aMCatNLO gridpack at NLO is naturally implementing multi-thread. LO gridpacks are, however, generated by a more aggregated mode—the "gridpack" mode (be aware that the name conflicts with CMS nomenclature), defined by the
gridpack
parameter. This mode initially supports single-thread only to generate LHE events. Since v2.6.1, a new feature is developed by the MG5_aMCatNLO team to support threaded event-production from the "gridpack" mode. The new
scriptName
, as specified above, provides an interface to the MG5_aMCatNLO threaded workflow.
The philosophy for this multi-thread implementation is as follows. When producing LHE events from an MG5_aMCatNLO "gridpack" generated by the "gridpack" mode, intermediate files are released to the "gridpack". Therefore, it naturally does not support multithreading due to the file overlapping across the threads.
In MG5_aMCatNLO v2.6.1, the "gridpack" framework is updated by the team to allow the "gridpack" to be stored in a read-only filesystem, while still enable the event generation. During the LHE production, the program detected the filesystem is read-only, hence writes the intermediate results elsewhere. Harnessing this feature, event generation can be launched from different external folders as multiple instances, while avoiding file overlapping. This is the main idea to develop the above script for multithreading.
Please note that the script can only cooperate with the CMS gridpack using MG5_aMCatNLO with version >= 2.6.1, in the LO mode. Errors will be reported if unqualified gridpack is detected.
Use case
This method can help save computing time for extremely large gridpack or complex ME process. However,
validation on the computing and physics performance is suggested before real production is launched. Some case studies on real samples are provided in the last section.
Pythia8ConcurrentHadronizerFilter
The
Pythia8HadronizerFilter
module hadronizes the EDM LHE events using
PYTHIA8
, then produces the EDM GEN events. A multi-thread–enhanced version is designed and employed, for this module that does not include a
ExternalDecays
parameter. To activate this feature, you can simply switch the module name to
Pythia8ConcurrentHadronizerFilter
. The configuration looks like
generator = cms.EDFilter("Pythia8ConcurrentHadronizerFilter",
...
)
Underlying concept
The original
Pythia8HadronizerFilter
is inefficient in event production, as only one event will process at a time. For example, in a case we specify 4 threads, the main instance assigns one event to a thread, while the other 3 threads have to wait until the event has finished processing. The main instance collects the event information, then assigns a second event to probably another thread, and the procedure goes on.
Pythia8ConcurrentHadronizerFilter
enables all threads to produce event concurrently. The main instance assigns new event to the thread once the previous event in this thread has finished. Click the button for more details.
The reason that
Pythia8HadronizerFilter
produce events in a less efficient way, is that the base class
HadronizerFilter
inherits from the
one::EDFilter
class. This prototype is one of the threaded framework documented in
MultithreadedFrameworkDesignDiscussions. It deals with one event at a time as we describe above to ensure the thread-safety. As some generator classes used in CMSSW are not thread-safe, which means, declaration of multiple instances of a generator class may cause conflicts between each other, the general
HadronizerFilter
is thus using the thread-safe but less efficient way to process events. On the other hands, Pythia8 class is a thread-safe class—multiple
Pythia
instances can process concurrently in different threads, therefore a dedicated base class
ConcurrentHadronizerFilter
is designed based on the
global::EDFilter
prototype that allows multiple threads running in a concurrent mode.
Pythia8ConcurrentHadronizerFilter
is the class based on
ConcurrentHadronizerFilter
that runs one Pythia8 class instance on each thread.
A
HadronizerFilter
module processes event in two steps. First, a "hadronizer" (Pythia8) help to hadronize the LHE to produce GEN events. Then, if a "decayer" is specified (e.g. EvtGen, Tauola), the GEN event is further processed by the "decayer". As "decayers" are not thread-safe, the
HadronizerFilter
takes a conservative step to process one event at a time, as mentioned above.
ConcurrentHadronizerFilter
runs the hadronizer at multiple threads, and does not include any "decayers".
Use case
Every
Pythia8HadronizerFilter
module that does not contain the
ExternalDecays
parameter is suggested to be replaced by
Pythia8ConcurrentHadronizerFilter
.
ExternalGeneratorFilter
A
GeneratorFilter
base module is similar to
HadronizerFilter
but produces events directly from a general-purpose generator, e.g. Pythia8, Herwig7, or Sherpa. For all modules based on
GeneratorFilter
, CMSSW provides the "wrapper"
ExternalGeneratorFilter
that enable multi-threaded event production. The configuration goes as follows (taking
Pythia8GeneratorFilter
as an example).
_generator = cms.EDFilter("Pythia8GeneratorFilter",
...
)
from GeneratorInterface.Core.ExternalGeneratorFilter import ExternalGeneratorFilter
generator = ExternalGeneratorFilter(_generator)
Note that you can either set the fragment or the python configuration file in this way. You may alternatively wrap the
cms.EDFilter
object by
ExternalGeneratorFilter
directly without the need of intermediate variable
_generator
. If you stick to the above way to set the fragment, be sure to start the intermediate variable name with an underline so
cmsDriver.py
knows to treat it as an intermediate variable.
You can also apply the setting by adding the following argument to
cmsDriver.py
, while using the original fragment.
--customise_commands='from GeneratorInterface.Core.ExternalGeneratorFilter import ExternalGeneratorFilter; process.generator = ExternalGeneratorFilter(process.generator)'
Underlying concept
The
ExternalGeneratorFilter
implements multi-thread by controlling external
cmsRun
instances, as implied from its name. When N threads are specified, each one of N
cmsRun
thread (act as a controller) starts an external
cmsRun
instance (as a worker) to process events. The controller transit the event information and random number properly to the worker thread, the later produce one event and transit the event product back to the controller. This way, the thread-safety issue can be avoided, since the controller and worker belong to separate programs, only with a proper share of memory.
The philosophy is to use boost c++ interprocess library to set off the interprocess communication between controller
cmsRun
and the worker
cmsRun
. The workers are created at the beginning of the program and will sleep unless they are called for single-event production. During the production, the controller rests and wait for the event to finish. Hence, the CPU resources are well-shared across all processes, and all threads can produce events simultaneously. The maximum per-event wait time for the controller is set to 60 s. It is doable to reset the time (e.g. to 120 s) by providing a parameter
_external_process_waitTime_=cms.untracked.uint32(120)
to the
ExternalGeneratorFilter
module.
Use case
All modules based on the
GeneratorFilter
class that generate events at the beginning stage (not using the previously generated LHE events) can use this method to enable multithreading.
Case study
In what follows we present four case studies for different use of the multithreading features.
Use ExternalGeneratorFilter with EvtGen
We first look at a physics process Λ
b0 → Λμμ events produced by EvtGen. The production uses the
Pythia8GeneratorFilter
module with
EvtGen
as the "decayer".
The process in the default settings suffers from low filter efficiency and cannot pass the McM validation. One solution is to use
ExternalGeneratorFilter
to launch the
Pythia8GeneratorFilter
concurrently to increase the event production speed. The fragment implementing this feature can be found here
LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_cff.py.
A full recipe that produces 1000 GEN events on 4 threads is shown as follows.
export SCRAM_ARCH=slc7_amd64_gcc700
source /cvmfs/cms.cern.ch/cmsset_default.sh
if [ -r CMSSW_10_6_19/src ] ; then
echo release CMSSW_10_6_19 already exists
else
scram p CMSSW CMSSW_10_6_19
fi
cd CMSSW_10_6_19/src
eval `scram runtime -sh`
curl -s -k https://twiki.cern.ch/twiki/pub/Sandbox/TestTopic11111266/LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_cff.py.txt \
--retry 3 --create-dirs -o Configuration/GenProduction/python/LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_cff.py
scram b
cd ../..
cmsDriver.py Configuration/GenProduction/python/LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_cff.py \
--python_filename LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_1_cfg.py \
--eventcontent RAWSIM --customise Configuration/DataProcessing/Utils.addMonitoring --datatier GEN \
--fileout file:LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen__GEN.root \
--conditions 106X_mc2017_realistic_v6 --beamspot Realistic25ns13TeVEarly2017Collision \
--customise_commands process.source.numberEventsInLuminosityBlock="cms.untracked.uint32(100)" \
--step GEN --geometry DB:Extended --era Run2_2017 --no_exec --mc -n 1000 --nThreads 4
cmsRun -e -j LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_report.xml LbToLambda0MuMu_Mufilter_SoftQCDnonD_TuneCP5_13TeV-pythia8-evtgen_1_cfg.py
Use ExternalGeneratorFilter with pure Pythia8
This example shows the Drell-Yan process with μμ final state produced by Pythia8 generator alone, using
Pythia8GeneratorFilter
. Similar to the above example, by wrapping the module with
ExternalGeneratorFilter
the production speed can be improved. The fragment implementing this feature can be found here
DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_cff.py.
A full recipe that produces 1000 GEN events on 4 threads is shown as follows.
export SCRAM_ARCH=slc7_amd64_gcc700
source /cvmfs/cms.cern.ch/cmsset_default.sh
if [ -r CMSSW_10_6_19/src ] ; then
echo release CMSSW_10_6_19 already exists
else
scram p CMSSW CMSSW_10_6_19
fi
cd CMSSW_10_6_19/src
eval `scram runtime -sh`
curl -s -k https://twiki.cern.ch/twiki/pub/Sandbox/TestTopic11111266/DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_cff.py.txt \
--retry 3 --create-dirs -o Configuration/GenProduction/python/DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_cff.py
scram b
cd ../..
cmsDriver.py Configuration/GenProduction/python/DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_cff.py \
--python_filename DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_1_cfg.py \
--eventcontent RAWSIM --customise Configuration/DataProcessing/Utils.addMonitoring --datatier GEN \
--fileout file:DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8__GEN.root \
--conditions 106X_mc2017_realistic_v6 --beamspot Realistic25ns13TeVEarly2017Collision \
--customise_commands process.source.numberEventsInLuminosityBlock="cms.untracked.uint32(100)" \
--step GEN --geometry DB:Extended --era Run2_2017 --no_exec --mc -n 1000 --nThreads 4
cmsRun -e -j DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_report.xml DYToMuMu_M-200To400_NNPDF30NLO_TuneCUETP8M1_13TeV-pythia8_1_cfg.py
Use concurrent ExternalLHEProducer for Powheg gridpack with Pythia8ConcurrentHadronizerFilter
This example shows the semileptonic tt̅ events produced by Powheg and hadronized by Pythia8. We consider using the concurrent mode of
ExternalLHEProducer
since the Powheg gridpack is pretty light-weighted. We also replace the Pythia8 hadronization module
Pythia8ConcurrentHadronizerFilter
as the module does not contain any decayers. The fragment implementing the two features are shown in
TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8_cff.py.
A full recipe that produces 1000 GEN events on 4 threads is shown as follows.
source /cvmfs/cms.cern.ch/cmsset_default.sh
if [ -r CMSSW_10_6_19/src ] ; then
echo release CMSSW_10_6_19 already exists
else
scram p CMSSW CMSSW_10_6_19
fi
cd CMSSW_10_6_19/src
eval `scram runtime -sh`
curl -s -k https://twiki.cern.ch/twiki/pub/Sandbox/TestTopic11111266/TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8_cff.py.txt \
--retry 3 --create-dirs -o Configuration/GenProduction/python/TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8_cff.py
scram b
cd ../..
cmsDriver.py Configuration/GenProduction/python/TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8_cff.py \
--python_filename TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8_1_cfg.py \
--eventcontent LHE,RAWSIM --customise Configuration/DataProcessing/Utils.addMonitoring --datatier LHE,GEN \
--fileout file:TTToSemiLeptonic_TuneCP5_13TeV-powheg-pythia8__wmLHEGEN.root \
--conditions 106X_mc2017_realistic_v6 --beamspot Realistic25ns13TeVEarly2017Collision \
--customise_commands process.source.numberEventsInLuminosityBlock="cms.untracked.uint32(100)" \
--step LHE,GEN --geometry DB:Extended --era Run2_2017 --no_exec --mc -n 1000 --nThreads 4
cmsRun -e -j YYY_report.xml YYY_1_cfg.py
Use concurrent ExternalLHEProducer for MG5_aMCatNLO NLO gridpack with Pythia8ConcurrentHadronizerFilter
This example shows gg → H+012j events produced by MG5_aMCatNLO at NLO mode with FxFx matching and hadronized by Pythia8. Given that the ME NLO calculation may be complicated and the RAM usage may become the bottleneck for multithreading, the concurrent ExternalLHEProducer may not outperform the default configuration, where MG5_aMCatNLO uses its own multi-threaded event generation strategy at NLO. Local tests are suggested to evaluate the two methods.
However, the Pythia8 hadronization module can be replaced to
Pythia8ConcurrentHadronizerFilter
as always recommended. The fragment implementing both concurrent
ExternalLHEProducer
and
Pythia8ConcurrentHadronizerFilter
are shown in
GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_cff.py.
A full recipe that produces 1000 GEN events on 4 threads is shown as follows.
export SCRAM_ARCH=slc7_amd64_gcc700
source /cvmfs/cms.cern.ch/cmsset_default.sh
if [ -r CMSSW_10_6_19/src ] ; then
echo release CMSSW_10_6_19 already exists
else
scram p CMSSW CMSSW_10_6_19
fi
cd CMSSW_10_6_19/src
eval `scram runtime -sh`
curl -s -k https://twiki.cern.ch/twiki/pub/Sandbox/TestTopic11111266/GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_cff.py.txt \
--retry 3 --create-dirs -o Configuration/GenProduction/python/GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_cff.py
scram b
cd ../..
cmsDriver.py Configuration/GenProduction/python/GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_cff.py \
--python_filename GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_1_cfg.py \
--eventcontent LHE,RAWSIM --customise Configuration/DataProcessing/Utils.addMonitoring --datatier LHE,GEN \
--fileout file:GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8__wmLHEGEN.root \
--conditions 106X_mc2017_realistic_v6 --beamspot Realistic25ns13TeVEarly2017Collision \
--customise_commands process.source.numberEventsInLuminosityBlock="cms.untracked.uint32(100)" \
--step LHE,GEN --geometry DB:Extended --era Run2_2017 --no_exec --mc -n 1000 --nThreads 4
cmsRun -e -j GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_report.xml GluGluHToGG_M126_TuneCP5_13TeV-amcatnloFXFX-pythia8_1_cfg.py
Use internal multithreading for MG5_aMCatNLO LO gridpack with Pythia8ConcurrentHadronizerFilter
This example shows W → ℓν+1234j events produced by MG5_aMCatNLO at LO with MLM matching and hadronized by Pythia8. Since this process has a very complex jet multiplicity, the gridpack contains complicated phase-space information and tends to be large. The large I/O consumption is the main concern during event generation. Therefore, we consider using the internal multithreading feature in MG5_aMCatNLO at LO mode. Again, local validations are suggested to compare the performance of two multithreading techniques in the LHE production.
Similarly, the Pythia8 hadronization module is replaced to
Pythia8ConcurrentHadronizerFilter
as recommended. The fragment implementing the internal solution to MG5_aMCatNLO at LO production and
Pythia8ConcurrentHadronizerFilter
is shown in
WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_cff.py.
A full recipe that produces 1000 GEN events on 4 threads is shown as follows.
export SCRAM_ARCH=slc7_amd64_gcc700
source /cvmfs/cms.cern.ch/cmsset_default.sh
if [ -r CMSSW_10_6_19/src ] ; then
echo release CMSSW_10_6_19 already exists
else
scram p CMSSW CMSSW_10_6_19
fi
cd CMSSW_10_6_19/src
eval `scram runtime -sh`
curl -s -k https://twiki.cern.ch/twiki/pub/Sandbox/TestTopic11111266/WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_cff.py.txt \
--retry 3 --create-dirs -o Configuration/GenProduction/python/WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_cff.py
scram b
cd ../..
cmsDriver.py Configuration/GenProduction/python/WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_cff.py \
--python_filename WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_1_cfg.py \
--eventcontent LHE,RAWSIM --customise Configuration/DataProcessing/Utils.addMonitoring --datatier LHE,GEN \
--fileout file:WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8__wmLHEGEN.root \
--conditions 106X_mc2017_realistic_v6 --beamspot Realistic25ns13TeVEarly2017Collision \
--customise_commands process.source.numberEventsInLuminosityBlock="cms.untracked.uint32(100)" \
--step LHE,GEN --geometry DB:Extended --era Run2_2017 --no_exec --mc -n 1000 --nThreads 4
cmsRun -e -j WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_report.xml WJetsToLNu_TuneCP5_13TeV-madgraph261MLM-pythia8_1_cfg.py