SpECTRE  v2025.03.17
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
RelativisticEuler::Solutions::RotatingStar Class Reference

A solution obtained by reading in rotating neutron star initial data from the RotNS code based on and . More...

#include <RotatingStar.hpp>

Classes

struct  PolytropicConstant
 The polytropic constant of the fluid. More...
 
struct  RotNsFilename
 The path to the RotNS data file. More...
 

Public Types

using options = implementation defined
 
- Public Types inherited from RelativisticEuler::AnalyticSolution< 3 >
using tags = implementation defined
 

Public Member Functions

 RotatingStar (const RotatingStar &)
 
RotatingStaroperator= (const RotatingStar &)
 
 RotatingStar (RotatingStar &&)=default
 
RotatingStaroperator= (RotatingStar &&)=default
 
 RotatingStar (std::string rot_ns_filename, std::unique_ptr< EquationsOfState::EquationOfState< true, 1 > > equation_of_state)
 
 RotatingStar (std::string rot_ns_filename, double polytropic_constant)
 
auto get_clone () const -> std::unique_ptr< evolution::initial_data::InitialData > override
 
template<typename DataType , typename... Tags>
tuples::TaggedTuple< Tags... > variables (const tnsr::I< DataType, 3 > &x, const double, tmpl::list< Tags... >) const
 Retrieve a collection of variables at (x, t)
 
void pup (PUP::er &p) override
 
const EquationsOfState::EquationOfState< true, 1 > & equation_of_state () const
 
double equatorial_radius () const
 
virtual auto get_clone () const -> std::unique_ptr< InitialData >=0
 
- Public Member Functions inherited from hydro::TemperatureInitialization< RotatingStar >
auto variables (const tnsr::I< DataType, Dim > &x, tmpl::list< hydro::Tags::Temperature< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::Temperature< DataType > >
 
auto variables (const tnsr::I< DataType, Dim > &x, const double t, tmpl::list< hydro::Tags::Temperature< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::Temperature< DataType > >
 
auto variables (ExtraVars &extra_variables, const tnsr::I< DataType, Dim > &x, Args &... extra_args, tmpl::list< hydro::Tags::Temperature< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::Temperature< DataType > >
 

Static Public Attributes

static constexpr Options::String help
 
- Static Public Attributes inherited from RelativisticEuler::AnalyticSolution< 3 >
static constexpr size_t volume_dim
 

Protected Types

template<typename DataType >
using DerivLapse = ::Tags::deriv< gr::Tags::Lapse< DataType >, tmpl::size_t< 3 >, Frame::Inertial >
 
template<typename DataType >
using DerivShift = ::Tags::deriv< gr::Tags::Shift< DataType, 3 >, tmpl::size_t< 3 >, Frame::Inertial >
 
template<typename DataType >
using DerivSpatialMetric = ::Tags::deriv< gr::Tags::SpatialMetric< DataType, 3 >, tmpl::size_t< 3 >, Frame::Inertial >
 

Protected Member Functions

template<typename DataType >
void interpolate_vars_if_necessary (gsl::not_null< IntermediateVariables< DataType > * > vars) const
 
template<typename DataType >
void interpolate_deriv_vars_if_necessary (gsl::not_null< IntermediateVariables< DataType > * > vars) const
 
template<typename DataType >
Scalar< DataType > lapse (const DataType &gamma, const DataType &rho) const
 
template<typename DataType >
tnsr::I< DataType, 3, Frame::Inertialshift (const DataType &omega, const DataType &phi, const DataType &radius, const DataType &sin_theta) const
 
template<typename DataType >
tnsr::ii< DataType, 3, Frame::Inertialspatial_metric (const DataType &gamma, const DataType &rho, const DataType &alpha, const DataType &phi) const
 
template<typename DataType >
tnsr::II< DataType, 3, Frame::Inertialinverse_spatial_metric (const DataType &gamma, const DataType &rho, const DataType &alpha, const DataType &phi) const
 
template<typename DataType >
Scalar< DataType > sqrt_det_spatial_metric (const DataType &gamma, const DataType &rho, const DataType &alpha) const
 
template<typename DataType >
auto make_metric_data (size_t num_points) const -> typename IntermediateVariables< DataType >::MetricData
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::RestMassDensity< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::RestMassDensity< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::ElectronFraction< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::ElectronFraction< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::SpecificEnthalpy< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::SpecificEnthalpy< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::Temperature< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::Temperature< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::Pressure< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::Pressure< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::SpecificInternalEnergy< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::SpecificInternalEnergy< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::SpatialVelocity< DataType, 3 > >) const -> tuples::TaggedTuple< hydro::Tags::SpatialVelocity< DataType, 3 > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::LorentzFactor< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::LorentzFactor< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::MagneticField< DataType, 3 > >) const -> tuples::TaggedTuple< hydro::Tags::MagneticField< DataType, 3 > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< hydro::Tags::DivergenceCleaningField< DataType > >) const -> tuples::TaggedTuple< hydro::Tags::DivergenceCleaningField< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::Lapse< DataType > >) const -> tuples::TaggedTuple< gr::Tags::Lapse< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::Shift< DataType, 3 > >) const -> tuples::TaggedTuple< gr::Tags::Shift< DataType, 3 > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::SpatialMetric< DataType, 3 > >) const -> tuples::TaggedTuple< gr::Tags::SpatialMetric< DataType, 3 > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::SqrtDetSpatialMetric< DataType > >) const -> tuples::TaggedTuple< gr::Tags::SqrtDetSpatialMetric< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::InverseSpatialMetric< DataType, 3 > >) const -> tuples::TaggedTuple< gr::Tags::InverseSpatialMetric< DataType, 3 > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list<::Tags::dt< gr::Tags::Lapse< DataType > > >) const -> tuples::TaggedTuple<::Tags::dt< gr::Tags::Lapse< DataType > > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< DerivLapse< DataType > >) const -> tuples::TaggedTuple< DerivLapse< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list<::Tags::dt< gr::Tags::Shift< DataType, 3 > > >) const -> tuples::TaggedTuple<::Tags::dt< gr::Tags::Shift< DataType, 3 > > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< DerivShift< DataType > >) const -> tuples::TaggedTuple< DerivShift< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list<::Tags::dt< gr::Tags::SpatialMetric< DataType, 3 > > >) const -> tuples::TaggedTuple<::Tags::dt< gr::Tags::SpatialMetric< DataType, 3 > > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< DerivSpatialMetric< DataType > >) const -> tuples::TaggedTuple< DerivSpatialMetric< DataType > >
 
template<typename DataType >
auto variables (gsl::not_null< IntermediateVariables< DataType > * > vars, const tnsr::I< DataType, 3 > &x, tmpl::list< gr::Tags::ExtrinsicCurvature< DataType, 3 > >) const -> tuples::TaggedTuple< gr::Tags::ExtrinsicCurvature< DataType, 3 > >
 

Protected Attributes

std::string rot_ns_filename_ {}
 
detail::CstSolution cst_solution_ {}
 
double polytropic_constant_ = std::numeric_limits<double>::signaling_NaN()
 
double polytropic_exponent_ = std::numeric_limits<double>::signaling_NaN()
 
bool is_polytrope_ {}
 
std::unique_ptr< EquationsOfState::EquationOfState< true, 1 > > equation_of_state_
 

Static Protected Attributes

static constexpr double atmosphere_floor_ = 1.e-50
 

Friends

bool operator== (const RotatingStar &lhs, const RotatingStar &rhs)
 

Detailed Description

A solution obtained by reading in rotating neutron star initial data from the RotNS code based on and .

The code that generates the initial data is part of a private SXS repository called RotNS.

The metric in spherical coordinates is given by

(1)ds2=eγ+ρdt2+e2α(dr2+r2dθ2)+eγρr2sin2(θ)(dϕωdt)2.

We use rotation about the z-axis. That is,

(2)gtt=eγ+ρ+eγρr2sin2(θ)ω2(3)grr=e2α(4)gθθ=e2αr2(5)gϕϕ=eγρr2sin2(θ)(6)gtϕ=eγρr2sin2(θ)ω.

We can transform from spherical to Cartesian coordinates using

(7)(r,θ,ϕ)(x,y,z)=(cos(ϕ)sin(θ)sin(θ)sin(ϕ)cos(θ)cos(ϕ)cos(θ)rsin(ϕ)cos(θ)rsin(θ)rsin(ϕ)rsin(θ)cos(ϕ)rsin(θ)0)

and

(8)(x,y,z)(r,θ,ϕ)=(cos(ϕ)sin(θ)rcos(ϕ)cos(θ)rsin(θ)sin(ϕ)sin(θ)sin(ϕ)rsin(ϕ)cos(θ)rcos(ϕ)sin(θ)cos(θ)rsin(θ)0)

We denote the lapse as N since α is already being used,

(9)N=e(γ+ρ)/2.

The shift is

(10)βϕ=ω

so

(11)βx=ϕxβϕ=sin(ϕ)ωrsin(θ),(12)βy=ϕyβϕ=cos(ϕ)ωrsin(θ),(13)βz=0.

The spatial metric is

γxx=sin2(θ)cos2(ϕ)e(2α)+cos2(θ)cos2(ϕ)e(2α)+sin2(ϕ)e(γρ)(14)=cos2(ϕ)e2α+sin2(ϕ)e(γρ)γxy=sin2(θ)cos(ϕ)sin(ϕ)e(2α)+cos2(θ)cos(ϕ)sin(ϕ)e(2α)sin(ϕ)cos(ϕ)e(γρ)(15)=cos(ϕ)sin(ϕ)e2αsin(ϕ)cos(ϕ)e(γρ)γyy=sin2(θ)sin2(ϕ)e(2α)+cos2(θ)sin2(ϕ)e(2α)+cos2(ϕ)e(γρ)(16)=sin2(ϕ)e2α+cos2(ϕ)e(γρ)(17)γxz=sin(θ)cos(ϕ)cos(θ)e2αcos(θ)cos(ϕ)sin(θ)e2α=0(18)γyz=sin(θ)sin(ϕ)cos(θ)e2αcos(θ)sin(ϕ)sin(θ)e2α=0(19)γzz=cos2(θ)e2α+sin2(θ)e2α=e2α

and its determinant is

(20)γ=e4α+(γρ)=e4αe(γρ).

At r=0 we have 2α=γρ and so the γxx=γyy=γzz=e2α and all other components are zero. The inverse spatial metric is given by

γxx=γyye2αe(γρ)=[sin2(ϕ)e2α+cos2(ϕ)e(γρ)]e2αe(γρ)(21)=sin2(ϕ)e(γρ)+cos2(ϕ)e2αγyy=γxxe2αe(γρ)=[cos2(ϕ)e2α+sin2(ϕ)e(γρ)]e2αe(γρ)(22)=cos2(ϕ)e(γρ)+sin2(ϕ)e2αγxy=γxye2αe(γρ)=[cos(ϕ)sin(ϕ)e2αsin(ϕ)cos(ϕ)e(γρ)]e2αe(γρ)=cos(ϕ)sin(ϕ)e(γρ)sin(ϕ)cos(ϕ)e2α(23)=cos(ϕ)sin(ϕ)[e2αe(γρ)](24)γxz=0(25)γyz=0(26)γzz=e2α.

The 4-velocity in spherical coordinates is given by

(27)ua¯=e(ρ+γ)/21v2[1,0,0,Ω],

where

(28)v=(Ωω)rsin(θ)eρ.

Transforming to Cartesian coordinates we have

(29)ut=e(ρ+γ)/21v2(30)ux=ϕxuϕ=rsin(θ)sin(ϕ)utΩ(31)uy=ϕyuϕ=rsin(θ)cos(ϕ)utΩ(32)uz=0.

The Lorentz factor is given by

W=Nut=e(γ+ρ)/2e(ρ+γ)/21v2(33)=11v2.

Using

(34)vi=1N(uiut+βi)

we get

(35)vx=e(γ+ρ)/2rsin(θ)sin(ϕ)(Ωω)=e(γρ)/2sin(ϕ)v(36)vy=e(γ+ρ)/2rsin(θ)cos(ϕ)(Ωω)=e(γρ)/2cos(ϕ)v(37)vz=0.

Lowering with the spatial metric we get

vx=γxxvx+γxyvy=[cos2(ϕ)e2α+sin2(ϕ)eγρ]e(γρ)/2sin(ϕ)v+[cos(ϕ)sin(ϕ)e2αsin(ϕ)cos(ϕ)eγρ]e(γρ)/2cos(ϕ)v=e(γρ)/2vsin(ϕ)eγρ[sin2(ϕ)+cos2(ϕ)](38)=e(γρ)/2vsin(ϕ)vy=γyxvx+γyyvy=[cos(ϕ)sin(ϕ)e2αsin(ϕ)cos(ϕ)e(γρ)]e(γρ)/2sin(ϕ)v+[sin2(ϕ)e2α+cos2(ϕ)e(γρ)]e(γρ)/2cos(ϕ)v(39)=e(γρ)/2vcos(ϕ)(40)vz=0.

This is consistent with the Lorentz factor read off from ut since vivi=v2. For completeness, ui=Wvi so

(41)ux=e(γρ)/2vsin(ϕ)1v2(42)uy=e(γρ)/2vcos(ϕ)1v2(43)uz=0.

Warning
Near (within 1e-2) r=0 the numerical errors from interpolation and computing the metric derivatives by finite difference no longer cancel out and so the tilde_s time derivative only vanishes to roughly 1e-8 rather than machine precision. Computing the Cartesian derivatives from analytic differentiation of the radial and angular polynomial fits might improve the situation but is a decent about of work to implement.

Member Function Documentation

◆ get_clone()

auto RelativisticEuler::Solutions::RotatingStar::get_clone ( ) const -> std::unique_ptr< evolution::initial_data::InitialData >
overridevirtual

Member Data Documentation

◆ help

constexpr Options::String RelativisticEuler::Solutions::RotatingStar::help
staticconstexpr
Initial value:
= {
"Rotating neutron star initial data solved by the RotNS solver. The data "
"is read in from disk. If RotNS used a polytropic equation of state, use "
"the PolytropicConstant option, so that the polytropic index may be read "
"directly from the RotNS output. Otherwise, set the equation of state to "
"whichever one you would like to use to load the initial data. Note that "
"if a polytrope is put into the EquationOfState option rather than using "
"the PolytropicConstant option, the generic unit conversion will be used "
"instead of the specific polytrope one."}

The documentation for this class was generated from the following file: