SpECTRE  v2025.03.17
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Particles::MonteCarlo Namespace Reference

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

Namespaces

namespace  Tags
 Items related to the evolution of particles Items related to Monte-Carlo radiation transport Tags for MC.
 

Classes

struct  CellLightCrossingTimeCompute
 
struct  DiffusionMonteCarloParameters
 Precomputed quantities useful for the diffusion approximation in high-scattering opacity regions. We follow [77] 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  InertialFrameEnergyDensityCompute
 
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  McGhostZoneData
 The container for ghost zone data to be communicated in the Monte-Carlo algorithm. More...
 
struct  McGhostZoneDataInboxTag
 Inbox tag to be used before MC step. More...
 
class  MonteCarloOptions
 
struct  MortarData
 Structure used to gather ghost zone data for Monte-Carlo evolution. We need the rest mass density, electron fraction, temperature, and an estimate of the light-crossing time one cell deep within each neighboring element. More...
 
class  NeutrinoInteractionTable
 Class responsible for reading neutrino-matter interaction tables. More...
 
struct  Packet
 Struct representing a single Monte Carlo packet of neutrinos. More...
 
struct  System
 
struct  TemplatedLocalFunctions
 Structure containing Monte-Carlo function templated on EnergyBins and/or NeutrinoSpecies. More...
 
struct  TimeStepMutator
 Mutator advancing neutrinos by a single step. More...
 

Enumerations

enum class  CommunicationStep { PreStep , PostStep }
 This class is used to template communication actions in the Monte-Carlo code, to know which of the two communication steps we need to take care of (just before or just after the MC step). More...
 

Functions

void cell_light_crossing_time (gsl::not_null< Scalar< DataVector > * > cell_light_crossing_time, const Mesh< 3 > &mesh, const tnsr::I< DataVector, 3, Frame::Inertial > &inertial_coordinates, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3, Frame::Inertial > &shift, const tnsr::II< DataVector, 3, Frame::Inertial > &inv_spatial_metric)
 
void cell_proper_four_volume_finite_difference (gsl::not_null< Scalar< DataVector > * > cell_proper_four_volume, const Scalar< DataVector > &lapse, const Scalar< DataVector > &sqrt_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, size_t extended_idx, 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)
 
void inertial_frame_energy_density (gsl::not_null< Scalar< DataVector > * > fluid_frame_energy_density, const std::vector< Packet > &packets, const Scalar< DataVector > &lapse, const Scalar< DataVector > &sqrt_determinant_spatial_metric, const Mesh< 3 > &mesh, const Scalar< DataVector > &det_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, size_t extended_idx, 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. More...
 

Detailed Description

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

Enumeration Type Documentation

◆ CommunicationStep

This class is used to template communication actions in the Monte-Carlo code, to know which of the two communication steps we need to take care of (just before or just after the MC step).

Enumerator
PreStep 

PreStep should occur just before the MC step. It should send to the ghost zones of each element the fluid and metric data.

PostStep 

PostStep should occur just after a MC step. It sends packets that have left their current element to the element they now belong to. NOT CODE YET: This is also when information about energy/momentum/lepton number exchance in the ghost zones should be communicated back to the corresponding live point for coupling to the fluid evolution.

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,
size_t  extended_idx,
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 [76], 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

(1)E=Wαptγijuipj

◆ diffuse_packet()

void Particles::MonteCarlo::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,
size_t  extended_idx,
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.

The coupling terms include ghost zone information, while other tensors only include live points. extended_idx is the index of the packet in DataVectors including ghost zones.