IAP GITLAB
Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
corsika
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Issue analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Air Shower Physics
corsika
Commits
99521d16
Commit
99521d16
authored
1 year ago
by
Felix Riehn
Committed by
Alan Coleman
1 year ago
Browse files
Options
Downloads
Patches
Plain Diff
printouts
parent
9fd8732d
No related branches found
No related tags found
1 merge request
!594
Neutrino interactions with pythia 8310
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
corsika/detail/modules/pythia8/NeutrinoInteraction.inl
+19
-12
19 additions, 12 deletions
corsika/detail/modules/pythia8/NeutrinoInteraction.inl
corsika/modules/pythia8/NeutrinoInteraction.hpp
+3
-1
3 additions, 1 deletion
corsika/modules/pythia8/NeutrinoInteraction.hpp
with
22 additions
and
13 deletions
corsika/detail/modules/pythia8/NeutrinoInteraction.inl
+
19
−
12
View file @
99521d16
...
...
@@ -20,9 +20,14 @@ namespace corsika::pythia8 {
,
handle_nc_
(
handleNC
)
,
handle_cc_
(
handleCC
)
,
pythiaMain_
{
CORSIKA_Pythia8_XML_DIR
,
false
}
{
Pythia8
::
RndmEngine
*
rndm
=
new
corsika
::
pythia8
::
Random
();
pythiaMain_
.
setRndmEnginePtr
(
rndm
);
CORSIKA_LOG_INFO
(
"configure Pythia for primary neutrino interactions. NC={}, CC={}"
,
handle_nc_
,
handle_cc_
);
CORSIKA_LOG_INFO
(
"minimal Q2 in DIS: {} GeV2"
,
minQ2_
/
1
_GeV
/
1
_GeV
);
// CORSIKA_LOG_INFO("Pythia8 MPI init file: {}", mpiInitFile.native());
// Main Pythia object for managing the cascade evolution.
// Can also do decays, but no hard processes.
...
...
@@ -105,27 +110,25 @@ namespace corsika::pythia8 {
FourMomentum
const
&
projectileP4
,
FourMomentum
const
&
targetP4
)
{
CORSIKA_LOG_INFO
(
"P
ythia::NeutrinoInteraction: doInteraction: {} E={}GeV"
,
projectileId
,
projectileP4
.
getTimeLikeComponent
()
/
1
_GeV
);
CORSIKA_LOG_INFO
(
"P
rimary {} - p interaction with E_nu = {} GeV"
,
projectileId
,
projectileP4
.
getTimeLikeComponent
()
/
1
_GeV
);
// allow just one call to NeutrinoInteraction
if
(
!
count_
)
{
double
eProton
=
Proton
::
mass
/
1
_GeV
;
// 0.93827;
double
eElectron
=
projectileP4
.
getTimeLikeComponent
()
/
1
_GeV
;
// 270.5;
double
Q2min
=
minQ2_
/
1
_GeV
;
double
const
eProton
=
Proton
::
mass
/
1
_GeV
;
// 0.93827;
double
const
eElectron
=
projectileP4
.
getTimeLikeComponent
()
/
1
_GeV
;
// 270.5;
double
const
Q2min
=
minQ2_
/
1
_GeV
/
1
_GeV
;
int
const
idNow
=
static_cast
<
int
>
(
get_PDG
(
projectileId
));
// Set up incoming beams, for frame with unequal beam energies.
// projectile is along -z
pythiaMain_
.
readString
(
"Beams:frameType = 2"
);
// BeamA = proton.
// BeamA = proton.
(+z)
pythiaMain_
.
readString
(
"Beams:idA = 2212"
);
pythiaMain_
.
settings
.
parm
(
"Beams:eA"
,
eProton
);
// BeamB = electron.
int
const
idNow
=
static_cast
<
int
>
(
get_PDG
(
projectileId
));
pythiaMain_
.
readString
(
"Beams:idB = 12"
);
// pythiaMain_.settings.parm("Beams:idA", idNow);
// BeamB = neutrino (-z direction)
pythiaMain_
.
settings
.
mode
(
"Beams:idB"
,
idNow
);
pythiaMain_
.
settings
.
parm
(
"Beams:eB"
,
eElectron
);
// Set up DIS process within some phase space.
...
...
@@ -175,7 +178,7 @@ namespace corsika::pythia8 {
MomentumVector
Plab_final
{
labFrameBoost
.
getOriginalCS
()};
auto
Elab_final
=
HEPEnergyType
::
zero
();
CORSIKA_LOG_INFO
(
"p
ythia stack size {}"
,
eventMain
.
size
()
);
CORSIKA_LOG_INFO
(
"p
articles generated in neutrino interaction:"
);
for
(
int
i
=
0
;
i
<
eventMain
.
size
();
++
i
)
{
if
(
eventMain
[
i
].
isFinal
())
{
auto
const
&
p8p
=
eventMain
[
i
];
...
...
@@ -195,6 +198,8 @@ namespace corsika::pythia8 {
auto
pnew
=
view
.
addSecondary
(
std
::
make_tuple
(
pyId
,
Ekin
,
pyPlab
.
normalized
()));
CORSIKA_LOG_INFO
(
"id = {}, E = {} GeV, p = {} GeV"
,
pyId
,
Ekin
/
1
_GeV
,
pyPlab
.
getComponents
()
/
1
_GeV
);
Plab_final
+=
pnew
.
getMomentum
();
Elab_final
+=
pnew
.
getEnergy
();
}
...
...
@@ -249,6 +254,8 @@ namespace corsika::pythia8 {
"Elab_final= {}"
", Plab_final= {}"
,
Elab_final
/
1
_GeV
,
(
Plab_final
/
1
_GeV
).
getComponents
());
}
else
{
CORSIKA_LOG_ERROR
(
"Only one primary neutrino interaction allowed!"
);
}
count_
++
;
}
...
...
This diff is collapsed.
Click to expand it.
corsika/modules/pythia8/NeutrinoInteraction.hpp
+
3
−
1
View file @
99521d16
...
...
@@ -20,6 +20,8 @@
namespace
corsika
::
pythia8
{
using
HEPEnergyTypeSqr
=
decltype
(
1
_GeV
*
1
_GeV
);
class
NeutrinoInteraction
:
public
InteractionProcess
<
NeutrinoInteraction
>
{
public:
...
...
@@ -68,7 +70,7 @@ namespace corsika::pythia8 {
bool
const
print_listing_
=
false
;
bool
const
handle_nc_
;
bool
const
handle_cc_
;
HEPEnergyType
const
minQ2_
=
25
_GeV
;
HEPEnergyType
Sqr
minQ2_
=
25
_GeV
*
1
_GeV
;
Pythia8
::
Pythia
pythiaMain_
;
};
}
// namespace corsika::pythia8
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment