SpECTRE  v2024.05.11
Particles::MonteCarlo Namespace Reference

Items related to the evolution of particles Items related to Monte-Carlo radiation transport. More...

Classes

struct  DiffusionMonteCarloParameters
 Precomputed quantities useful for the diffusion approximation in high-scattering opacity regions. We follow [70] Note that r_d in that manuscript should be (distance diffused)/(time elapsed) i.e. the paper is missing a normalization of r_d by delta_t. The upper bound of the integral in Eq (30) should just be sqrt(3*tau/4) All quantities in this struct are fixed; we could also just hard-code them to save us the on-the-fly calculation, or have a single instance of this struct in the global cache. More...
 
struct  InverseJacobianInertialToFluidCompute
 Inverse Jacobian of the map from inertial coordinate to an orthonormal frame comoving with the fluid. That frame uses the 4-velocity as its time axis, and constructs the other members of the tetrads using Gram-Schmidt's algorithm. More...
 
struct  Packet
 Struct representing a single Monte Carlo packet of neutrinos. More...
 
struct  TemplatedLocalFunctions
 Structure containing Monte-Carlo function templated on EnergyBins and/or NeutrinoSpecies. More...
 

Functions

void cell_proper_four_volume_finite_difference (gsl::not_null< Scalar< DataVector > * > cell_proper_four_volume, const Scalar< DataVector > &lapse, const Scalar< DataVector > &determinant_spatial_metric, double time_step, const Mesh< 3 > &mesh, const Scalar< DataVector > &det_jacobian_logical_to_inertial)
 Proper 4-volume of a cell for a time step time_step. We assume that determinant_spatial_metric is given in inertial coordinates, hence the need for det_jacobian_logical_to_inertial.
 
void cell_inertial_coordinate_three_volume_finite_difference (gsl::not_null< Scalar< DataVector > * > cell_inertial_three_volume, const Mesh< 3 > &mesh, const Scalar< DataVector > &det_jacobian_logical_to_inertial)
 3-volume of a cell in inertial coordinate. Note that this is the coordinate volume, not the proper volume. This quantity is needed as a normalization factor for the terms coupling Monte-Carlo transport to the fluid evolution (as we evolved densitized fluid variables in inertial coordinates).
 
void AddCouplingTermsForPropagation (gsl::not_null< Scalar< DataVector > * > coupling_tilde_tau, gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * > coupling_tilde_s, gsl::not_null< Scalar< DataVector > * > coupling_rho_ye, const Packet &packet, double dt, double absorption_opacity, double scattering_opacity, double fluid_frame_energy, double lapse, double lorentz_factor, const std::array< double, 3 > &lower_spatial_four_velocity_packet)
 Contribution to the neutrino-matter coupling terms from evolving a packet by time dt, given absorption and scattering opacities. More...
 
void evolve_single_packet_on_geodesic (gsl::not_null< Packet * > packet, double time_step, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3, Frame::Inertial > &shift, const tnsr::i< DataVector, 3, Frame::Inertial > &d_lapse, const tnsr::iJ< DataVector, 3, Frame::Inertial > &d_shift, const tnsr::iJJ< DataVector, 3, Frame::Inertial > &d_inv_spatial_metric, const tnsr::II< DataVector, 3, Frame::Inertial > &inv_spatial_metric, const std::optional< tnsr::I< DataVector, 3, Frame::Inertial > > &mesh_velocity, const InverseJacobian< DataVector, 3, Frame::ElementLogical, Frame::Inertial > &inverse_jacobian_logical_to_inertial)
 
double compute_fluid_frame_energy (const Packet &packet, const Scalar< DataVector > &lorentz_factor, const tnsr::i< DataVector, 3, Frame::Inertial > &lower_spatial_four_velocity, const Scalar< DataVector > &lapse, const tnsr::II< DataVector, 3, Frame::Inertial > &inv_spatial_metric)
 
void DiffusionPrecomputeForElement (gsl::not_null< DataVector * > prefactor_diffusion_time_vector, gsl::not_null< DataVector * > prefactor_diffusion_four_velocity, gsl::not_null< DataVector * > prefactor_diffusion_time_step, const Scalar< DataVector > &lorentz_factor, const tnsr::i< DataVector, 3, Frame::Inertial > &lower_spatial_four_velocity, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3, Frame::Inertial > &shift, const tnsr::ii< DataVector, 3, Frame::Inertial > &spatial_metric)
 Precompute quantities needed for evolving packets in the diffusion approximation, i.e. the transport velocity in logical coordinates and the coefficients defined by Eq (31-33)
 
std::array< double, 4 > scatter_packet (gsl::not_null< Packet * > packet, gsl::not_null< std::mt19937 * > random_number_generator, const double &fluid_frame_energy, const Jacobian< DataVector, 4, Frame::Inertial, Frame::Fluid > &inertial_to_fluid_jacobian, const InverseJacobian< DataVector, 4, Frame::Inertial, Frame::Fluid > &inertial_to_fluid_inverse_jacobian)
 Elastic scattering, which is just redrawing the momentum from an isotropic distribution in the fluid frame, at constant fluid frame energy.
 
void diffuse_packet (gsl::not_null< Packet * > packet, gsl::not_null< std::mt19937 * > random_number_generator, gsl::not_null< double * > neutrino_energy, gsl::not_null< Scalar< DataVector > * > coupling_tilde_tau, gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * > coupling_tilde_s, gsl::not_null< Scalar< DataVector > * > coupling_rho_ye, double time_step, const DiffusionMonteCarloParameters &diffusion_params, double absorption_opacity, double scattering_opacity, const Scalar< DataVector > &lorentz_factor, const tnsr::i< DataVector, 3, Frame::Inertial > &lower_spatial_four_velocity, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3, Frame::Inertial > &shift, const tnsr::i< DataVector, 3, Frame::Inertial > &d_lapse, const tnsr::iJ< DataVector, 3, Frame::Inertial > &d_shift, const tnsr::iJJ< DataVector, 3, Frame::Inertial > &d_inv_spatial_metric, const tnsr::ii< DataVector, 3, Frame::Inertial > &spatial_metric, const tnsr::II< DataVector, 3, Frame::Inertial > &inv_spatial_metric, const std::optional< tnsr::I< DataVector, 3, Frame::Inertial > > &mesh_velocity, const InverseJacobian< DataVector, 3, Frame::ElementLogical, Frame::Inertial > &inverse_jacobian_logical_to_inertial, const Jacobian< DataVector, 4, Frame::Inertial, Frame::Fluid > &inertial_to_fluid_jacobian, const InverseJacobian< DataVector, 4, Frame::Inertial, Frame::Fluid > &inertial_to_fluid_inverse_jacobian, const DataVector &prefactor_diffusion_time_step, const DataVector &prefactor_diffusion_four_velocity, const DataVector &prefactor_diffusion_time_vector)
 Evolve a packet for dt = time_step, assuming that we can use the diffusion approximation.
 

Detailed Description

Items related to the evolution of particles Items related to Monte-Carlo radiation transport.

Function Documentation

◆ AddCouplingTermsForPropagation()

void Particles::MonteCarlo::AddCouplingTermsForPropagation ( gsl::not_null< Scalar< DataVector > * >  coupling_tilde_tau,
gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * >  coupling_tilde_s,
gsl::not_null< Scalar< DataVector > * >  coupling_rho_ye,
const Packet packet,
double  dt,
double  absorption_opacity,
double  scattering_opacity,
double  fluid_frame_energy,
double  lapse,
double  lorentz_factor,
const std::array< double, 3 > &  lower_spatial_four_velocity_packet 
)

Contribution to the neutrino-matter coupling terms from evolving a packet by time dt, given absorption and scattering opacities.

Details

We calculate total momentum exchanges, not densities; thus when coupling to the evolution of the energy/momentum density, division by the spatial volume is necessary. We calculate the source terms of Eqs (62-64) of [69], with integrals calculated as in Eqs (6-10) of that manuscript (except that we do not divide by the coordinate volume V).

◆ compute_fluid_frame_energy()

double Particles::MonteCarlo::compute_fluid_frame_energy ( const Packet packet,
const Scalar< DataVector > &  lorentz_factor,
const tnsr::i< DataVector, 3, Frame::Inertial > &  lower_spatial_four_velocity,
const Scalar< DataVector > &  lapse,
const tnsr::II< DataVector, 3, Frame::Inertial > &  inv_spatial_metric 
)

Calculate energy of neutrinos in a frame comoving with the fluid

\begin{align} E = W \alpha p^t - \gamma^{ij} u_i p_j \end{align}