-
Notifications
You must be signed in to change notification settings - Fork 180
Add generator for LambdaB to deuteron #2237
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Open
ercolessi
wants to merge
3
commits into
AliceO2Group:master
Choose a base branch
from
ercolessi:master
base: master
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
+140
−0
Open
Changes from all commits
Commits
Show all changes
3 commits
Select commit
Hold shift + click to select a range
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,7 @@ | ||
| #NEV_TEST> 10 | ||
| ### The external generator derives from GeneratorPythia8. | ||
| [GeneratorExternal] | ||
| fileName=${O2DPG_MC_CONFIG_ROOT}/MC/config/PWGHF/external/generator/generator_pythia8_hfhadron_to_nuclei.C | ||
| funcName=generateHFHadToNuclei(3, {5122}, {1000010020}, true, 0.5) | ||
| [GeneratorPythia8] | ||
| config=${O2DPG_MC_CONFIG_ROOT}/MC/config/PWGHF/pythia8/generator/pythia8_lambdab_probQQtoQ.cfg | ||
57 changes: 57 additions & 0 deletions
57
MC/config/PWGHF/ini/tests/GeneratorHFTrigger_LambdaBToDeuteron.C
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,57 @@ | ||
| int External() | ||
| { | ||
| std::string path{"o2sim_Kine.root"}; | ||
| std::vector<int> checkPdgHadron{5122}; // Lambda_b | ||
| std::vector<int> nucleiDauPdg{1000010020}; // d, 3He, 3H | ||
|
|
||
| TFile file(path.c_str(), "READ"); | ||
| if (file.IsZombie()) | ||
| { | ||
| std::cerr << "Cannot open ROOT file " << path << "\n"; | ||
| return 1; | ||
| } | ||
|
|
||
| auto tree = (TTree *)file.Get("o2sim"); | ||
| std::vector<o2::MCTrack> *tracks{}; | ||
| tree->SetBranchAddress("MCTrack", &tracks); | ||
|
|
||
| int nSignals{}, nSignalGoodDecay{}; | ||
| auto nEvents = tree->GetEntries(); | ||
|
|
||
| for (int i = 0; i < nEvents; i++) | ||
| { | ||
| tree->GetEntry(i); | ||
| for (auto &track : *tracks) | ||
| { | ||
| auto pdg = track.GetPdgCode(); | ||
| if (std::find(checkPdgHadron.begin(), checkPdgHadron.end(), std::abs(pdg)) != checkPdgHadron.end()) // found signal | ||
| { | ||
| // count signal PDG | ||
| if(std::abs(track.GetRapidity()) > 1.5) continue; // skip if outside rapidity window | ||
| nSignals++; | ||
| for (int j{track.getFirstDaughterTrackId()}; j <= track.getLastDaughterTrackId(); ++j) | ||
| { | ||
| auto pdgDau = tracks->at(j).GetPdgCode(); | ||
| if (std::find(nucleiDauPdg.begin(), nucleiDauPdg.end(), std::abs(pdgDau)) != nucleiDauPdg.end()) | ||
| { | ||
| nSignalGoodDecay++; | ||
| } | ||
| } | ||
| } | ||
| } | ||
| } | ||
| std::cout << "--------------------------------\n"; | ||
| std::cout << "# Events: " << nEvents << "\n"; | ||
| std::cout <<"# signal hadrons: " << nSignals << "\n"; | ||
| std::cout <<"# signal hadrons decaying into nuclei: " << nSignalGoodDecay << "\n"; | ||
|
|
||
| float fracForcedDecays = nSignals ? float(nSignalGoodDecay) / nSignals : 0.0f; | ||
| float uncFracForcedDecays = nSignals ? std::sqrt(fracForcedDecays * (1 - fracForcedDecays) / nSignals) : 1.0f; | ||
| if (1 - fracForcedDecays > 0.2 + uncFracForcedDecays) // we put some tolerance (lambdaB in MB events do not coalesce) | ||
| { | ||
| std::cerr << "Fraction of signals decaying into nuclei: " << fracForcedDecays << ", lower than expected\n"; | ||
| return 1; | ||
| } | ||
|
|
||
| return 0; | ||
| } |
76 changes: 76 additions & 0 deletions
76
MC/config/PWGHF/pythia8/generator/pythia8_lambdab_probQQtoQ.cfg
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,76 @@ | ||
| ### beams | ||
| Beams:idA 2212 # proton | ||
| Beams:idB 2212 # proton | ||
| Beams:eCM 13600. # GeV | ||
|
|
||
| ### processes | ||
| SoftQCD:inelastic on # all inelastic processes | ||
|
|
||
| ### decays | ||
| ParticleDecays:limitTau0 on | ||
| ParticleDecays:tau0Max 10. | ||
|
|
||
| # enable deuteron production by coalescence collisions | ||
| HadronLevel:DeuteronProduction = on | ||
| DeuteronProduction:channels = {2212 2112 > 22} | ||
| DeuteronProduction:models = {0} | ||
| DeuteronProduction:norm = 1 | ||
| DeuteronProduction:parms = {0.5 1} # coalescence momentum p0 in GeV/c | ||
|
|
||
| ### switching on Pythia Mode2 | ||
| ColourReconnection:mode 1 | ||
| ColourReconnection:allowDoubleJunRem off | ||
| ColourReconnection:m0 0.3 | ||
| ColourReconnection:allowJunctions on | ||
| ColourReconnection:junctionCorrection 1.20 | ||
| ColourReconnection:timeDilationMode 2 | ||
| ColourReconnection:timeDilationPar 0.18 | ||
| StringPT:sigma 0.335 | ||
| StringZ:aLund 0.36 | ||
| StringZ:bLund 0.56 | ||
| #StringFlav:probQQtoQ 0.078 | ||
| StringFlav:probQQtoQ 0.17 | ||
| StringFlav:ProbStoUD 0.2 | ||
| StringFlav:probQQ1toQQ0join 0.0275,0.0275,0.0275,0.0275 | ||
| MultiPartonInteractions:pT0Ref 2.15 | ||
| BeamRemnants:remnantMode 1 | ||
| BeamRemnants:saturation 5 | ||
|
|
||
| # Correct decay lengths (wrong in PYTHIA8 decay table) | ||
| # Lb | ||
| 5122:tau0 = 0.4390 | ||
| # Xic0 | ||
| 4132:tau0 = 0.0455 | ||
| # OmegaC | ||
| 4332:tau0 = 0.0803 | ||
|
|
||
| ## Lb decay modes | ||
| 5122:onMode = off | ||
| 5122:oneChannel = 1 9.03e-05 0 2212 2212 -2212 -2212 2112 | ||
| 5122:addChannel = 1 0.00181 0 2212 2212 -2212 -2212 2112 111 | ||
| 5122:addChannel = 1 0.0181 0 2212 2212 -2212 -2212 2112 111 111 | ||
| 5122:addChannel = 1 0.0271 0 2212 2212 -2212 -2212 2112 211 -211 | ||
| 5122:addChannel = 1 0.0181 0 2212 2212 -2212 -2212 2112 211 211 -211 -211 | ||
| 5122:addChannel = 1 0.181 0 2212 2212 -2212 -2212 2112 211 -211 111 | ||
| 5122:addChannel = 1 0.199 0 2212 2212 -2212 -2212 2112 211 -211 111 111 | ||
| 5122:addChannel = 1 0.0271 0 2212 2212 -2212 -2212 2112 211 211 -211 -211 111 | ||
| 5122:addChannel = 1 0.00452 0 2212 2212 -2212 -2212 2112 211 111 | ||
| 5122:addChannel = 1 0.0361 0 2212 2212 2112 -2212 -2112 -211 111 | ||
| 5122:addChannel = 1 0.0452 0 2212 2212 2112 -2212 -2112 -211 111 111 | ||
| 5122:addChannel = 1 0.0542 0 2212 2212 2112 -2212 -2112 -211 -211 211 | ||
| 5122:addChannel = 1 0.361 0 2212 2212 2112 -2212 -2112 -211 -211 211 111 | ||
| 5122:addChannel = 1 0.0271 0 2212 2212 2112 -2212 -2112 -211 -211 211 111 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 111 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 211 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 211 211 211 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 211 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 211 111 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 211 211 211 111 | ||
| 5122:onIfMatch = 2212 2212 2212 2212 2112 211 111 | ||
| 5122:onIfMatch = 2212 2212 2112 2212 2112 211 111 | ||
| 5122:onIfMatch = 2212 2212 2112 2212 2112 211 111 111 | ||
| 5122:onIfMatch = 2212 2212 2112 2212 2112 211 211 211 | ||
| 5122:onIfMatch = 2212 2212 2112 2212 2112 211 211 211 111 | ||
| 5122:onIfMatch = 2212 2212 2112 2212 2112 211 211 211 111 111 |
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
Uh oh!
There was an error while loading. Please reload this page.