Skip to content

Adding new alternate neutrino cooling #1777

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
wants to merge 16 commits into
base: development
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions .github/workflows/good_defines.txt
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@ MICROPHYSICS_DEBUG
NAUX_NET
NETWORK_SOLVER
NEUTRINOS
NEUTRINO_METHOD
NEW_NETWORK_IMPLEMENTATION
NONAKA_PLOT
NSE
Expand Down
2 changes: 1 addition & 1 deletion .github/workflows/test_neutrinos.yml
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ jobs:
run: |
cd unit_test/test_neutrino_cooling
./main3d.gnu.ex inputs amrex.fpe_trap_{invalid,zero,overflow}=1
../../external/amrex/Tools/Plotfile/fextrema.gnu.ex test_sneut5 > test.out
../../external/amrex/Tools/Plotfile/fextrema.gnu.ex test_neutrino.sneut5 > test.out

- name: Print backtrace
if: ${{ failure() && hashFiles('unit_test/test_neutrino_cooling/Backtrace.0') != '' }}
Expand Down
9 changes: 9 additions & 0 deletions Make.Microphysics_extern
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,15 @@ else
$(error Invalid value for SCREEN_METHOD)
endif

NEUTRINO_METHOD ?= kipp
ifeq ($(NEUTRINO_METHOD), kipp)
DEFINES += -DNEUTRINO_METHOD=NEUTRINO_METHOD_kipp
else ifeq ($(NEUTRINO_METHOD), sneut5)
DEFINES += -DNEUTRINO_METHOD=NEUTRINO_METHOD_sneut5
else
$(error Invalid value for NEUTRINO_METHOD: $(NEUTRINO_METHOD))
endif

ifeq ($(EOS_DIR), gamma_law_general)
override EOS_DIR := gamma_law
$(warning gamma_law_general has been renamed gamma_law)
Expand Down
3 changes: 2 additions & 1 deletion neutrinos/Make.package
Original file line number Diff line number Diff line change
@@ -1 +1,2 @@
CEXE_headers += sneut5.H
CEXE_headers += neutrino.H
CEXE_headers += sneut5.H kipp.H
172 changes: 172 additions & 0 deletions neutrinos/kipp.H
Original file line number Diff line number Diff line change
@@ -0,0 +1,172 @@
// This routine calculates the neutrino cooling rates using the approximations from Kippenhahn et al. 1989
// using autodiff

#ifndef KIP_NEUT_H
#define KIP_NEUT_H

#include<iostream>
#include<AMReX_REAL.H>
#include<fundamental_constants.H>
#include<microphysics_autodiff.H>

using namespace amrex::literals;

template<typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
// pair annihilation
number_t kip_pair(const number_t& temp, const amrex::Real& rho,
const number_t& abar, const number_t& zbar){
number_t T9 = temp * 1.0e-9_rt;
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

add:

    amrex::ignore_unused(abar);
    amrex::ignore_unused(zbar);

to suppress unused var warnings

// eq. 18.81 - fixed using woosley lecture notes
if (T9 < 2.0e0_rt){
return (4.9e18_rt/rho) * amrex::Math::powi<3>(T9) * admath::exp(-11.86/T9);
}
else{
return (3.2e15_rt/rho) * amrex::Math::powi<9>(T9);
}
}

template<typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
// photoneutrinos
number_t kip_phot(const number_t& temp, const amrex::Real& rho,
const number_t& abar, const number_t& zbar){
number_t T9 = temp * 1.0e-9_rt;
// eq. 18.82
number_t ep_1 = (1.103e13_rt/rho) * amrex::Math::powi<9>(T9) * admath::exp(-5.93e0_rt / T9);
number_t ep_2 = 0.976e8_rt * amrex::Math::powi<8>(T9) / (1.0e0_rt + 4.2e0_rt * T9);
number_t rho_bar = 6.446e-6_rt * rho / (T9 + 4.2e0_rt * amrex::Math::powi<2>(T9));
number_t mu_e = abar / zbar;

return ep_1 + ep_2 / (mu_e + rho_bar);
}

template<typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
// plasmaneutrinos
number_t kip_plas(const number_t& temp, const amrex::Real& rho,
const number_t& abar, const number_t& zbar){
number_t n_e = rho / (C::m_u * abar/zbar);
// assuming non-degeneracy
// 18.83 replaced by aprox in woosley lecture notes
number_t w0 = admath::sqrt((4.0e0_rt * M_PI * amrex::Math::powi<2>(C::q_e) * n_e) / C::m_e);
number_t gamma = (C::hbar * w0) / (C::k_B * temp);
number_t lambda = (C::k_B * temp) / (C::m_e * amrex::Math::powi<2>(C::c_light));
// 18.85
if (gamma <= 1.0e0_rt) {
// return (7.4e21_rt * amrex::Math::powi<6>(C::hbar * w0/(C::m_e * amrex::Math::powi<2>(C::c_light)))
// * amrex::Math::powi<3>(lambda))/rho;
number_t epsilon_p = (C::hbar * w0) / (C::m_e * C::c_light * C::c_light);
return (7.4e21_rt * amrex::Math::powi<6>(epsilon_p) * amrex::Math::powi<3>(lambda)) / rho;

} else{
return (3.3e21_rt * admath::pow((C::hbar * w0/(C::m_e * amrex::Math::powi<2>(C::c_light))),7.5e0_rt)
* admath::pow(lambda, 3.0e0_rt/2.0e0_rt) * admath::exp(-gamma))/rho;
}

}

template<typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
// bremsstrahlung neutrinos
number_t kip_brem(const number_t& temp, const amrex::Real& rho,
const number_t& abar, const number_t& zbar){
// eq. 18.86
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

add:

    amrex::ignore_unused(rho);

number_t T8 = temp * 1.0e-8_rt;
return 0.76e0_rt * (amrex::Math::powi<2>(zbar)/abar) * amrex::Math::powi<6>(T8);
}

template<typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
number_t kipp(const number_t& temp, const amrex::Real rho,
const number_t& abar, const number_t& zbar, amrex::Real& pair,
amrex::Real& phot, amrex::Real& plas, amrex::Real& brem){

if (temp < 1.0e7_rt){
return 0.0e0_rt;
}

// calculating individual contributions of all processes
number_t kip_pair_dual = kip_pair(temp, rho, abar, zbar);
number_t kip_phot_dual = kip_phot(temp, rho, abar, zbar);
number_t kip_plas_dual = kip_plas(temp, rho, abar, zbar);
number_t kip_brem_dual = kip_brem(temp, rho, abar, zbar);

// extracting contributions to store in plotfile
pair = autodiff::val(kip_pair_dual);
phot = autodiff::val(kip_phot_dual);
plas = autodiff::val(kip_plas_dual);
brem = autodiff::val(kip_brem_dual);

// total neutrino cooling rate
number_t snu = kip_pair_dual + kip_phot_dual + kip_plas_dual + kip_brem_dual;

return snu;

}

template <int do_derivatives>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
void kipp(const amrex::Real& temp, const amrex::Real& rho, const amrex::Real& abar,
const amrex::Real& zbar, amrex::Real& snu, amrex::Real& dsnudt,
amrex::Real& dsnudrho, amrex::Real& dsnudz, amrex::Real& dsnuda,
amrex::Real& pair, amrex::Real& phot, amrex::Real& plas, amrex::Real& brem) {

/*
input:
temp = temperature
rho = density
abar = mean atomic weight
zbar = mean charge

output:
snu = total neutrino cooling rate in erg/g/s
dsnudt = derivative of snu with respect to temp
dsnudz = derivative of snu with respect to zbar
dsnuda = derivative of snu with respect to abar
pair = pair annihilation contribution in erg/g/s
phot = photoneutrino contribution in erg/g/s
plas = plasma neutrino contribution in erg/g/s
brem = bremsstrahlung contribution rate in erg/g/s
*/

// autodiff wrapper
if constexpr (do_derivatives){

using dual_t = autodiff::dual_array<1,3>;
dual_t temp_dual = temp;
dual_t abar_dual = abar;
dual_t zbar_dual = zbar;

autodiff::seed_array(temp_dual, abar_dual, zbar_dual);
dual_t snu_dual = kipp(temp_dual, rho, abar_dual, zbar_dual, pair, phot, plas, brem);

snu = autodiff::val(snu_dual);
const auto& grad = autodiff::derivative(snu_dual);
dsnudt = grad(1);
dsnudrho = 0.0e0_rt;
dsnuda = grad(2);
dsnudz = grad(3);
} else {
snu = kipp(temp, rho, abar, zbar, pair, phot, plas, brem);
dsnudt = 0.0e0_rt;
dsnudrho = 0.0e0_rt;
dsnuda = 0.0e0_rt;
dsnudz = 0.0e0_rt;
}
}

//overloading to pass arguments in all cases
template <int do_derivatives>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
void kipp(const amrex::Real temp, const amrex::Real den,
const amrex::Real abar, const amrex::Real zbar,
amrex::Real& snu, amrex::Real& dsnudt, amrex::Real& dsnudd,
amrex::Real& dsnuda, amrex::Real& dsnudz){

amrex::Real pair, phot, plas, brem;
kipp<do_derivatives>(temp, den, abar, zbar, snu, dsnudt, dsnudd, dsnuda,
dsnudz, pair, phot, plas, brem);
}

#endif // KIP_NEUT_H
39 changes: 39 additions & 0 deletions neutrinos/neutrino.H
Original file line number Diff line number Diff line change
@@ -0,0 +1,39 @@
#ifndef NEUTRINO_H
#define NEUTRINO_H

// Define available neutrino cooling methods
#define NEUTRINO_METHOD_kipp 1
#define NEUTRINO_METHOD_sneut5 2

#include <AMReX_REAL.H>
#include <AMReX_Array.H>
#include <microphysics_autodiff.H>
#include <fundamental_constants.H>

#include <kipp.H>
#include <sneut5.H>

using namespace amrex::literals;

#if NEUTRINO_METHOD == NEUTRINO_METHOD_kipp
const std::string neutrino_name = "kipp";
#elif NEUTRINO_METHOD == NEUTRINO_METHOD_sneut5
const std::string neutrino_name = "sneut5";
#endif

template <typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
number_t neutrino_cooling(const number_t& temp, const amrex::Real dens,
const number_t& abar, const number_t& zbar){

#if NEUTRINO_METHOD == NEUTRINO_METHOD_kipp
amrex::Real pair, phot, plas, brem;
return kipp(temp, dens, abar, zbar, pair, phot, plas, brem);

#elif NEUTRINO_METHOD == NEUTRINO_METHOD_sneut5
amrex::Real pair, phot, plas, brem;
return sneut5(temp, dens, abar, zbar, pair, phot, plas, brem);
#endif
}

#endif // NEUTRINO_H
28 changes: 24 additions & 4 deletions neutrinos/sneut5.H
Original file line number Diff line number Diff line change
Expand Up @@ -832,7 +832,8 @@ number_t nu_recomb(const sneutf_t<number_t>& sf) {
template <typename number_t>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
number_t sneut5(const number_t& temp, const amrex::Real den,
const number_t& abar, const number_t& zbar)
const number_t& abar, const number_t& zbar, amrex::Real& pair,
amrex::Real& phot, amrex::Real& plas, amrex::Real& brem)
{
/*
this routine computes thermal neutrino losses from the analytic fits of
Expand Down Expand Up @@ -872,6 +873,11 @@ number_t sneut5(const number_t& temp, const amrex::Real den,
sbrem *= sf.deni;
sreco *= sf.deni;

pair = autodiff::val(spair);
phot = autodiff::val(sphot);
plas = autodiff::val(splas);
brem = autodiff::val(sbrem);

// the total neutrino loss rate
number_t snu = splas + spair + sphot + sbrem;
if (neutrino_cooling_rp::include_recomb) {
Expand All @@ -885,7 +891,8 @@ AMREX_GPU_HOST_DEVICE AMREX_INLINE
void sneut5(const amrex::Real temp, const amrex::Real den,
const amrex::Real abar, const amrex::Real zbar,
amrex::Real& snu, amrex::Real& dsnudt, amrex::Real& dsnudd,
amrex::Real& dsnuda, amrex::Real& dsnudz)
amrex::Real& dsnuda, amrex::Real& dsnudz, amrex::Real& pair,
amrex::Real& phot, amrex::Real& plas, amrex::Real& brem)
{
// autodiff wrapper
if constexpr (do_derivatives) {
Expand All @@ -894,15 +901,16 @@ void sneut5(const amrex::Real temp, const amrex::Real den,
dual_t abar_dual = abar;
dual_t zbar_dual = zbar;
autodiff::seed_array(temp_dual, abar_dual, zbar_dual);
dual_t snu_dual = sneut5(temp_dual, den, abar_dual, zbar_dual);
dual_t snu_dual = sneut5(temp_dual, den, abar_dual, zbar_dual,
pair, phot, plas, brem);
snu = autodiff::val(snu_dual);
const auto& grad = autodiff::derivative(snu_dual);
dsnudt = grad(1);
dsnudd = 0.0e0_rt;
dsnuda = grad(2);
dsnudz = grad(3);
} else {
snu = sneut5(temp, den, abar, zbar);
snu = sneut5(temp, den, abar, zbar, pair, phot, plas, brem);
dsnudt = 0.0e0_rt;
dsnudd = 0.0e0_rt;
dsnuda = 0.0e0_rt;
Expand All @@ -911,4 +919,16 @@ void sneut5(const amrex::Real temp, const amrex::Real den,

}

template <int do_derivatives>
AMREX_GPU_HOST_DEVICE AMREX_INLINE
void sneut5(const amrex::Real temp, const amrex::Real den,
const amrex::Real abar, const amrex::Real zbar,
amrex::Real& snu, amrex::Real& dsnudt, amrex::Real& dsnudd,
amrex::Real& dsnuda, amrex::Real& dsnudz){

amrex::Real pair, phot, plas, brem;
sneut5<do_derivatives>(temp, den, abar, zbar, snu, dsnudt, dsnudd, dsnuda,
dsnudz, pair, phot, plas, brem);
}

#endif
3 changes: 3 additions & 0 deletions unit_test/test_neutrino_cooling/GNUmakefile
Original file line number Diff line number Diff line change
Expand Up @@ -26,6 +26,9 @@ EOS_DIR := helmholtz
# This sets the network directory
NETWORK_DIR := aprox21

# This sets the neurtino cooling method
NEUTRINO_METHOD := sneut5

# This isn't actually used but we need VODE to compile with CUDA
INTEGRATOR_DIR := VODE

Expand Down
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
plotfile = test_sneut5
plotfile = test_neutrino.sneut5
time = 0
variables minimum value maximum value
density 10 5000000000
Expand Down Expand Up @@ -28,4 +28,8 @@
dsneutdt 0 2.9114557633e+14
dsneutda -2.0028760264e+16 1.6034657416e+20
dsneutdz -1.832532276e+20 2.2890025857e+16
pair 0 3.3206297949e+23
phot 0 2.5493486735e+14
plas 0 7.284558679e+12
brem -2.5119883805e+12 83290903991

5 changes: 4 additions & 1 deletion unit_test/test_neutrino_cooling/main.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,9 @@ using namespace amrex;

#include <network.H>
#include <eos.H>

#include <neutrino.H>
#include <kipp.H>
#include <sneut5.H>

#include <variables.H>
Expand Down Expand Up @@ -194,7 +197,7 @@ void main_main ()
ParallelDescriptor::ReduceRealMax(stop_time, IOProc);


std::string name = "test_sneut5";
std::string name = "test_neutrino." + neutrino_name;

// Write a plotfile
WriteSingleLevelPlotfile(name, state, names, geom, time, 0);
Expand Down
Loading
Loading