sph_particles Module

This module contains the definition of TYPES, PROCEDURES and variables needed to set up the \(\mathrm{SPH}\) \(\mathrm{ID}\)



Uses

  • module~~sph_particles~~UsesGraph module~sph_particles sph_particles module~id_base id_base module~sph_particles->module~id_base module~utility utility module~sph_particles->module~utility timing timing module~sph_particles->timing module~id_base->module~utility module~id_base->timing constants constants module~utility->constants matrix matrix module~utility->matrix

Used by

  • module~~sph_particles~~UsedByGraph module~sph_particles sph_particles module~access access module~access->module~sph_particles module~apm apm module~apm->module~sph_particles module~bssn_formulation bssn_formulation module~bssn_formulation->module~sph_particles module~standard_tpo_formulation standard_tpo_formulation module~bssn_formulation->module~standard_tpo_formulation module~compose compose module~compose->module~sph_particles module~constructor_bin constructor_bin module~constructor_bin->module~sph_particles module~constructor_std constructor_std module~constructor_std->module~sph_particles module~ellipsoidal_surfaces ellipsoidal_surfaces module~ellipsoidal_surfaces->module~sph_particles module~handle_positions handle_positions module~handle_positions->module~sph_particles module~io io module~io->module~sph_particles module~lattices lattices module~lattices->module~sph_particles module~memory memory module~memory->module~sph_particles module~quality_indicators quality_indicators module~quality_indicators->module~sph_particles module~recovery recovery module~recovery->module~sph_particles module~redistribute_nu redistribute_nu module~redistribute_nu->module~sph_particles module~sph_variables sph_variables module~sph_variables->module~sph_particles module~standard_tpo_formulation->module~sph_particles program~convergence_test convergence_test program~convergence_test->module~sph_particles program~convergence_test->module~bssn_formulation program~convergence_test->module~standard_tpo_formulation module~cauchy_convergence_test cauchy_convergence_test program~convergence_test->module~cauchy_convergence_test program~sphincs_id sphincs_id program~sphincs_id->module~sph_particles program~sphincs_id->module~bssn_formulation module~access~9 access module~access~9->module~standard_tpo_formulation module~analysis analysis module~analysis->module~standard_tpo_formulation module~bssn_variables bssn_variables module~bssn_variables->module~bssn_formulation module~cauchy_convergence_test->module~standard_tpo_formulation module~constraints constraints module~constraints->module~bssn_formulation module~constructor constructor module~constructor->module~bssn_formulation module~io~2 io module~io~2->module~standard_tpo_formulation module~io~3 io module~io~3->module~bssn_formulation module~landau_lifshitz landau_lifshitz module~landau_lifshitz->module~bssn_formulation module~memory~2 memory module~memory~2->module~bssn_formulation module~recovery_m2p recovery_m2p module~recovery_m2p->module~standard_tpo_formulation module~ricci ricci module~ricci->module~bssn_formulation module~sph_adm_variables sph_adm_variables module~sph_adm_variables->module~standard_tpo_formulation module~standard_tpo_variables standard_tpo_variables module~standard_tpo_variables->module~standard_tpo_formulation module~perform_test perform_test module~perform_test->module~cauchy_convergence_test module~shared_grid shared_grid module~shared_grid->module~cauchy_convergence_test

Contents


Variables

Type Visibility Attributes Name Initial
integer, public, parameter :: eos$id = 1

First component of array eos_parameters for a polytropic \(\mathrm{EOS}\)

integer, public, parameter :: id_particles_from_formatted_file = 0

Identifier for a particle distribution read from formatted file

integer, public, parameter :: id_particles_on_ellipsoidal_surfaces = 2

Identifier for particle distribution on ellipsoidal surfaces

integer, public, parameter :: id_particles_on_lattice = 1

Identifier for a particle distribution on a lattice

integer, public, parameter :: max_it_tree = 1

When computing the neighbours' tree with the SUBROUTINE exact_nei_tree_update, it can happen that some particles do not have exactly 300 neighours. If this happens, the smoothing lenghts of the particles in the tree-cell are increased by a factor of 3, and used as a new guess to recompute the entire tree. max_it_tree specifies how many times the process should be iterated, if needed.

Note that, at the end of the iteration, the smoothing lengths are checked for them being NaNs or 0. For the particles at which this happens, the smoothing lengths are computed brute-force using the SUBROUTINE find_h_backup, so that such particles have exactly 300 neighbours.

integer, public, parameter :: poly$gamma = 2

Second component of array eos_parameters for a polytropic \(\mathrm{EOS}\)

integer, public, parameter :: poly$kappa = 3

Third component of array eos_parameters for a polytropic \(\mathrm{EOS}\)

integer, public, parameter, DIMENSION(4) :: pwp$gamma = [3, 4, 5, 6]

Third to sixth component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter, DIMENSION(4) :: pwp$kappa = [7, 8, 9, 10]

Seventh to tenth component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter :: pwp$log10p1 = 11

Eleventh component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter :: pwp$log10rho0 = 12

Twelfth component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter :: pwp$log10rho1 = 13

Thirteenth component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter :: pwp$log10rho2 = 14

Fourteenth component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)

integer, public, parameter :: pwp$npoly = 2

Second component of array eos_parameters for a piecewise polytropic \(\mathrm{EOS}\)


Interfaces

interface

  • public module subroutine allocate_particles_memory(this)

    Allocates allocatable arrays member of a particles object

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

interface

  • public module subroutine analyze_hydro(this, namefile)

    Scans the hydro fields taken from to look for negative or zero values

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    character(len=*), intent(inout), optional :: namefile

    Name of the formatted file where the particle positions at which some of the hydro fields are negative or zero are printed to

interface

  • public module function check_particle_position(npart, pos, pos_a) result(cnt)

    Return the number of times that pos_a appears in the array pos

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart
    double precision, intent(in), DIMENSION(3,npart) :: pos
    double precision, intent(in), DIMENSION(3) :: pos_a

    Return Value integer

interface

  • public module subroutine check_particle_positions(npart, pos, debug)

    Check that the particles are not at the same positions

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart

    Number of particles

    double precision, intent(in), DIMENSION(3,npart) :: pos

    Array of particle positions

    logical, intent(in), optional :: debug

    TRUE to debug the SUBROUTINE, FALSE otherwise

interface

  • public module subroutine compute_Ye(this)

    Interpolates linearly the electron fraction at the particle densities; that is, assigns at the particle positions

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

interface

  • public module subroutine compute_and_print_quality_indicators(npart, pos, h, nu, nstar, path)

    Compute the quality indicators, referring to

    Daniel J. Price, Smoothed Particle Hydrodynamics and Magnetohydrodynamics. Journal of Computational Physics, 231, 3, 759-794 (2012) DOI: 10.1016/j.jcp.2010.12.011 http://arxiv.org/abs/1012.1885 eqs.(64), (67) and (74-75)

    Rosswog, S. SPH Methods in the Modelling of Compact Objects. Living Rev Comput Astrophys 1, 1 (2015). https://doi.org/10.1007/lrca-2015-1. https://arxiv.org/abs/1406.4224. eqs.(6) and (9)

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart
    double precision, intent(in), DIMENSION(3,npart) :: pos
    double precision, intent(in), DIMENSION(npart) :: h
    double precision, intent(in), DIMENSION(npart) :: nu
    double precision, intent(in), DIMENSION(npart) :: nstar
    character(len=*), intent(in), optional :: path

    Path to which saving the output file containing the quality indicators.

interface

  • public module subroutine compute_and_print_sph_variables(this, namefile)

    Computes the SPH variables at the particle positions, and optionally prints them to a binary file to be read by and , and to a formatted file to be read by , by calling print_formatted_id_particles

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    character(len=*), intent(inout), optional :: namefile

    Name of the formatted file where the SPH ID is printed to

interface

  • public module subroutine compute_sph_hydro(this, npart_in, npart_fin, eqos, nlrf, u, Pr, enthalpy, cs, verbose)

    Computes the hydro fields on a section of the particles specified as input. First, computes the \(\mathrm{SPH}\) pressure starting from the \(\mathrm{SPH}\) baryon mass density, and the specific internal energy. The pressure is computed differently for different \(\mathrm{EOS}\), and for cold and hot systems. Then computes the enthalpy and the sound speed accordingly.

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: npart_in

    First index of the desired section of the particles

    integer, intent(in) :: npart_fin

    Last index of the desired section of the particles

    class(eos), intent(in) :: eqos

    \(\mathrm{EOS}\) to be used

    double precision, intent(in), DIMENSION(npart_fin - npart_in + 1) :: nlrf

    Baryon mass density in the local rest frame

    double precision, intent(inout), DIMENSION(npart_fin - npart_in + 1) :: u

    Specific internal energy

    double precision, intent(inout), DIMENSION(npart_fin - npart_in + 1) :: Pr

    Pressure

    double precision, intent(inout), DIMENSION(npart_fin - npart_in + 1) :: enthalpy

    Enthalpy

    double precision, intent(inout), DIMENSION(npart_fin - npart_in + 1) :: cs

    Speed of sound

    logical, intent(in), optional :: verbose

interface

  • public module function construct_particles_bin(id, namefile) result(parts)

    Constructs a particles object from an \(\mathrm{ID}\) binary file

    Arguments

    Type IntentOptional Attributes Name
    class(idbase), intent(inout) :: id

    idbase object representing the BNS for which we want to place particles

    character(len=*), intent(inout) :: namefile

    Name of the \(\mathrm{ID}\) binary file

    Return Value type(particles)

    Constructed particles object

interface

  • public module function construct_particles_std(id, dist) result(parts)

    Standard constructor for a particles object

    Arguments

    Type IntentOptional Attributes Name
    class(idbase), intent(inout) :: id

    idbase object representing the BNS for which we want to place particles

    integer, intent(in) :: dist

    Identifier of the desired particle distribution:

    • 0: Read particle positions (and optionally the baryon number per particle ) from a formatted file

    • 1: Place particles on several lattices, each one surrounding a matter object

    • 3: Place particles on ellipsoidal surfaces inside each matter object

    Return Value type(particles)

    Constructed particles object

interface

  • public module subroutine correct_center_of_mass(npart, pos, nu, get_density, validate_pos, com_star, verbose)

    Translate the particles so that their center of mass coincides with the center of mass of the star, given by \(\mathrm{ID}\)

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart
    double precision, intent(inout), DIMENSION(3,npart) :: pos
    double precision, intent(inout), DIMENSION(npart) :: nu
    function get_density(x, y, z) result(density)
    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x
    double precision, intent(in) :: y
    double precision, intent(in) :: z
    Return Value double precision
    function validate_pos(x, y, z) result(answer)
    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x
    double precision, intent(in) :: y
    double precision, intent(in) :: z
    Return Value logical
    double precision, intent(in), DIMENSION(3) :: com_star
    logical, intent(in), optional :: verbose

interface

  • public module subroutine deallocate_particles_memory(this)

    Deallocates allocatable arrays member of a particles object

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

interface

  • public module subroutine destruct_particles(this)

    Arguments

    Type IntentOptional Attributes Name
    type(particles), intent(inout) :: this

    Finalizer (Destructor) of particles object particles object which this PROCEDURE is a member of

interface

  • public module subroutine find_particles_above_xy_plane(npart, pos, npart_above_xy, above_xy_plane_a)

    Mirror the particle with z>0 with respect to the xy plane, to impose the equatorial-plane symmetry

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart
    double precision, intent(in), DIMENSION(3,npart) :: pos
    integer, intent(out) :: npart_above_xy
    integer, intent(out), DIMENSION(:), ALLOCATABLE :: above_xy_plane_a

interface

interface

  • public pure module function get_eos_id(this, i_matter) result(eos_id)

    Returns compose_eos

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: i_matter

    Index of the matter object

    Return Value integer

    (compose_eos)

interface

  • public pure module function get_g3(this) result(g3)

    Returns the spatial metric on the particles

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(6,this% npart)

    (g_xx,g_xy,g_xz,g_yy,g_yz,g_zz)

interface

  • public pure module function get_h(this) result(h)

    Returns h

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    h

interface

  • public pure module function get_lapse(this) result(lapse)

    Returns lapse

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(this% npart)

    lapse

interface

  • public pure module function get_n_matter(this) result(n_matter)

    Returns n_matter

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value integer

    n_matter

interface

  • public module subroutine get_neighbours_bf(ipart, npart, pos, h, dimensions, nnei, neilist)

    Get neighbours of particle ipart in a "brute force" way

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: ipart
    integer, intent(in) :: npart
    double precision, intent(in) :: pos(dimensions,npart)
    double precision, intent(in) :: h(npart)
    integer, intent(in) :: dimensions
    integer, intent(out) :: nnei
    integer, intent(out) :: neilist(npart)

interface

  • public pure module function get_nlrf(this) result(nlrf)

    Returns nlrf

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    nlrf

interface

  • public pure module function get_nlrf_sph(this) result(nlrf_sph)

    Returns nlrf_sph

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    nlrf_sph

interface

  • public pure module function get_npart(this) result(n_part)

    Returns npart

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value integer

    npart

interface

  • public pure module function get_npart_i(this, i_matter) result(n_part)

    Returns the number of particles on the object i_matter

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: i_matter

    Index of the matter object

    Return Value integer

    Number of particles on the object i_matter

interface

  • public pure module function get_nstar(this) result(nstar)

    Returns nstar

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    nstar

interface

  • public pure module function get_nstar_sph(this) result(nstar_sph)

    Returns nstar_sph

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    nstar_sph

interface

  • public pure module function get_nu(this) result(nu)

    Returns nu

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    nu

interface

  • public pure module function get_nuratio(this) result(nuratio)

    Returns nuratio

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision

    nuratio

interface

  • public pure module function get_nuratio_i(this, i_matter) result(nuratio)

    Returns the baryon number ratio on the object i_matter

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: i_matter

    Index of the matter object

    Return Value double precision

    Baryon number ratio on the object i_matter

interface

  • public pure module function get_object_of_particle(this, a)

    Returns the number of the matter object asociated with the particle number given as input. Example: give number as input; this particle number corresponds to a particle on matter object . This functions returns .

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: a

    Particle number

    Return Value integer

interface

  • public pure module function get_pos(this) result(pos_u)

    Returns pos

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:,:), ALLOCATABLE

    pos

interface

  • public pure module function get_pressure(this) result(pressure)

    Returns pressure

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    pressure

interface

  • public pure module function get_pressure_sph(this) result(pressure_sph)

    Returns pressure_sph

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    pressure_sph

interface

  • public pure module function get_shift(this) result(shift)

    Returns

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(3,this% npart)

    (shift_x,shift_y,shift_z)

interface

  • public pure module function get_theta(this) result(theta)

    Returns Theta

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    Theta

interface

interface

  • public pure module function get_u_sph(this) result(u_sph)

    Returns u_sph

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:), ALLOCATABLE

    u_sph

interface

  • public pure module function get_vel(this) result(vel)

    Returns v

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value double precision, DIMENSION(:,:), ALLOCATABLE

    v

interface

  • public module subroutine impose_equatorial_plane_symmetry(npart, pos, nu, com_star, verbose)

    Mirror the particle with z>0 with respect to the xy plane, to impose the equatorial-plane symmetry

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(inout) :: npart
    double precision, intent(inout), DIMENSION(3,npart) :: pos
    double precision, intent(inout), optional, DIMENSION(npart) :: nu
    double precision, intent(in), optional :: com_star
    logical, intent(in), optional :: verbose

interface

  • public module function is_empty(this) result(answer)

    Returns .TRUE if the particles object is empty, .FALSE otherwise

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(in) :: this

    particles object which this PROCEDURE is a member of

    Return Value logical

    .TRUE if the particles object is empty, .FALSE otherwise

public interface particles

Interface of TYPE particles

interface

  • public module subroutine perform_apm(get_density, get_nstar_id, get_pressure_id, compute_pressure, npart_output, pos_input, pvol, h_output, nu_output, center, com_star, mass, sizes, eqos, apm_max_it, max_inc, mass_it, correct_nu, nuratio_thres, nuratio_des, use_pressure, adapt_ghosts, move_away_ghosts, nx_gh, ny_gh, nz_gh, ghost_dist, use_atmosphere, remove_atmosphere, print_step, namefile_pos_id, namefile_pos, namefile_results, validate_position, surf)

    Performs the Artificial Pressure Method (APM)

    particles object which this PROCEDURE is a member of

    Arguments

    Type IntentOptional Attributes Name
    function get_density(x, y, z) result(density)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    Return Value double precision

    Baryon mass density at

    subroutine get_nstar_id(npart, x, y, z, nstar_sph, nstar_id, nlrf_sph, sqg)

    Computes the proper baryon number density at the particle positions

    Arguments
    Type IntentOptional Attributes Name
    integer, intent(in) :: npart

    Number of real particles (i.e., no ghost particles included here)

    double precision, intent(in) :: x(npart)

    Array of coordinates

    double precision, intent(in) :: y(npart)

    Array of coordinates

    double precision, intent(in) :: z(npart)

    Array of coordinates

    double precision, intent(in) :: nstar_sph(npart)

    \(\mathrm{SPH}\) proper baryon density

    double precision, intent(out) :: nstar_id(npart)

    Array to store the computed proper baryon number density

    double precision, intent(out) :: nlrf_sph(npart)

    Array to store the local rest frame baryon density computed from the \(\mathrm{SPH}\) proper baryon density

    double precision, intent(out) :: sqg(npart)

    Square root of minus the determinant of the spacetime metric

    function get_pressure_id(x, y, z) result(pressure)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    Return Value double precision

    Baryon mass density at

    subroutine compute_pressure(npart, x, y, z, nlrf, eqos, pressure, verbose)
    Arguments
    Type IntentOptional Attributes Name
    integer, intent(in) :: npart

    Returns the baryon mass density at the desired point

    double precision, intent(in) :: x(npart)

    coordinate of the desired point

    double precision, intent(in) :: y(npart)

    coordinate of the desired point

    double precision, intent(in) :: z(npart)

    coordinate of the desired point

    double precision, intent(in) :: nlrf(npart)

    Baryon mass density in the local rest frame

    type(eos), intent(in) :: eqos

    \(\mathrm{EOS}\) to use

    double precision, intent(inout) :: pressure(npart)

    Baryon mass density at

    logical, intent(in), optional :: verbose

    If .TRUE., print informative standard output about how the pressure is computed. Default is .TRUE.

    integer, intent(inout) :: npart_output

    Initial particle number

    double precision, intent(inout), DIMENSION(:,:), ALLOCATABLE :: pos_input

    Initial particle positions

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: pvol

    Initial particle volume

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: h_output

    Array to store the smoothing lengths computed at the end of the APM iteration

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: nu_output

    Array to store the baryon number per particle computed at the end of the APM iteration

    double precision, intent(in), DIMENSION(3) :: center

    Center of the matter object, from the \(\mathrm{ID}\)

    double precision, intent(inout), DIMENSION(3) :: com_star

    Center of mass of the matter object, from the \(\mathrm{ID}\)

    double precision, intent(in) :: mass

    Mass of the matter object

    double precision, intent(in), DIMENSION(6) :: sizes

    Sizes of the matter object

    type(eos), intent(in) :: eqos

    \(\mathrm{EOS}\) to use when computing the pressure

    integer, intent(in) :: apm_max_it

    Maximum number of APM iterations, irrespective of the EXIT condition

    integer, intent(in) :: max_inc

    Sets the EXIT condition: If the average over all the particles of the relative error in the density estimate grows max_inc times, exit the iteration.

    logical, intent(in) :: mass_it

    If .TRUE. performs a second iteration after the APM one, without moving the particles, changing their mass in order to better match the star density. The mass ratio grows very fast in all the tried experiments, hence the suggested value is .FALSE.

    logical, intent(in) :: correct_nu

    If .TRUE., the baryon number per particle nu is corrected to include the total baryonic masses of the stars.

    double precision, intent(in) :: nuratio_thres

    Maximum mass ratio (equivalently baryon number ratio) to be used in the one-time-only final correction of the particle masses to match the star density even better (without moving the particles)

    double precision, intent(in) :: nuratio_des

    Sets the EXIT condition: If the baryon number ratio is within 2.5% of nuratio_des, exit the iteration Set nuratio_des to 0 to deactivate and exit the APM iteration using max_inc

    logical, intent(in) :: use_pressure

    If .TRUE., uses the physical pressure computed with the \(\mathrm{EOS}\) using the SPH estimate of the density nlrf_sph, to compute the artificial pressure. Otherwise, the density variable nstar_sph is used to compute the artificial pressure

    logical, intent(in) :: adapt_ghosts

    If .TRUE., the ghost particles will be placed and have a baryon number such to reproduce the density of the outermost layers (r > 99% of the minimum radius) of the object. If .TRUE., the arguments nx_gh, ny_gh, nz_gh, ghost_dist are ignored; if .FALSE., they are instead used to place the ghost particles

    logical, intent(in) :: move_away_ghosts

    If .TRUE., the ghost particles will slowly move away from the surface of the star during the iteration (depending on the EXIT condition chosen, this happens in slightly different ways), to allow for the real particles to get closer to the surface

    integer, intent(in) :: nx_gh

    Number of lattice points in the x direction for ghosts

    integer, intent(in) :: ny_gh

    Number of lattice points in the y direction for ghosts

    integer, intent(in) :: nz_gh

    Number of lattice points in the z direction for ghosts

    double precision, intent(inout) :: ghost_dist

    Distance between the ghost particles and the surface of the matter object considered (star, ejecta, etc...)

    logical, intent(inout) :: use_atmosphere

    If .TRUE., allows particles to move where the density is 0, and displace them using only the smoothing term. This can be useful when the system has an irregular geometry, as, for example, an ejecta; .FALSE. otherwise

    logical, intent(inout) :: remove_atmosphere

    If .TRUE., removes the particles placed where the density is 0, at the end of the APM iteration; .FALSE. otherwise

    integer, intent(inout) :: print_step

    Prints the particle positions to a formatted file every print_step steps

    character(len=*), intent(inout), optional :: namefile_pos_id

    Name for the formatted file where the initial particle positions and the ghost positions will be printed

    character(len=*), intent(inout), optional :: namefile_pos

    Name for the formatted file where the particle positions and the ghost positions will be printed every 15 iterations

    character(len=*), intent(inout), optional :: namefile_results

    Name for the formatted file where various quantities related to the particle distribution, the baryon number particle and the kernel estimate of the density will be printed at the end of the APM iteration

    procedure(validate_position_int), optional :: validate_position

    Returns 1 if the position is not valid, 0 otherwise

    type(surface), intent(in), optional :: surf

    Surface of the matter object

interface

  • public module subroutine place_particles_ellipsoidal_surfaces(this, mass_star, radius, center, central_density, npart_des, npart_out, pos, pvol, nu, h, last_r, upper_bound, lower_bound, upper_factor, lower_factor, max_steps, filename_mass_profile, filename_shells_radii, filename_shells_pos, get_density, integrate_density, get_id, validate_position, radii, surf)

    Places particles on ellipsoidal surfaces on one star

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    double precision, intent(in) :: mass_star

    Baryonic mass of the star

    double precision, intent(in) :: radius

    Radius of the star in the x direction towards the companion

    double precision, intent(in), DIMENSION(3) :: center

    (x|) coordinate of the center of the star, i.e., of the point with highest density

    double precision, intent(in) :: central_density

    Central density of the star, i.e., highest density

    integer, intent(in) :: npart_des

    idbase object needed to access the BNS data

    integer, intent(out) :: npart_out

    Final number of particles on the star

    double precision, intent(inout), DIMENSION(:,:), ALLOCATABLE :: pos

    Array string the final positions

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: pvol

    Array soring the inal particle volumes

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: nu

    Array storing the final particle masses Array storing the particle baryon masses

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: h

    Array storing the initial guess for the particle smoothing lengths

    double precision, intent(in) :: last_r

    Radius of the last ellipsoidal surface

    double precision, intent(inout) :: upper_bound

    Desired upper bound for the differences between particle masses on neighbouring ellipsoidal surfaces

    double precision, intent(inout) :: lower_bound

    Desired lower bound for the differences between particle masses on neighbouring ellipsoidal surfaces

    double precision, intent(in) :: upper_factor

    If, after max_steps, the iteration did not converge, multiply upper_bound by upper_factor, and lower_bound by lower_factor. upper_factor >= 1, usually an increase of 1% works

    double precision, intent(in) :: lower_factor

    If, after max_steps, the iteration did not converge, multiply upper_bound by upper_factor, and lower_bound by lower_factor. lower_factor <= 1, usually a decrease of 1% works

    integer, intent(in) :: max_steps

    If, after max_steps, the iteration did not converge, multiply upper_bound by upper_factor, and lower_bound by lower_factor. max_steps >= 10. 100 is a nice value

    character(len=*), intent(inout), optional :: filename_mass_profile
    character(len=*), intent(inout), optional :: filename_shells_radii

    Name of the file to store the surface radii

    character(len=*), intent(inout), optional :: filename_shells_pos

    Name of the file to store the final particle positions

    function get_density(x, y, z) result(density)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    Return Value double precision

    Baryon mass density at

    subroutine integrate_density(center, radius, central_density, dr, dth, dphi, mass, mass_profile, mass_profile_idx, radii, surf)
    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in), DIMENSION(3) :: center

    Center of the star

    double precision, intent(in) :: radius

    Radius of the star

    double precision, intent(in) :: central_density

    Central density of the star

    double precision, intent(in) :: dr

    Integration steps

    double precision, intent(in) :: dth

    Integration steps

    double precision, intent(in) :: dphi

    Integration steps

    double precision, intent(inout) :: mass

    Integrated mass of the star

    double precision, intent(out), DIMENSION(3,0:NINT(radius/dr)) :: mass_profile

    Array storing the radial mass profile of the star

    integer, intent(out), DIMENSION(0:NINT(radius/dr)) :: mass_profile_idx

    Array to store the indices for array mass_profile, sorted so that mass_profile[mass_profile_idx] is in increasing order INTEGER, DIMENSION(:), ALLOCATABLE, INTENT(INOUT)::mass_profile_idx

    double precision, intent(in), optional, DIMENSION(2) :: radii
    type(surface), intent(in), optional :: surf

    Surface of the matter object

    subroutine get_id(x, y, z, sqdetg, baryon_density, gamma_euler)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    double precision, intent(out) :: sqdetg
    double precision, intent(out) :: baryon_density
    double precision, intent(out) :: gamma_euler
    procedure(validate_position_int), optional :: validate_position

    Returns 1 if the position is not valid, 0 otherwise

    double precision, intent(in), optional, DIMENSION(2) :: radii
    type(surface), intent(in), optional :: surf

    Surface of the matter object

interface

  • public module subroutine place_particles_lattice(this, central_density, xmin, xmax, ymin, ymax, zmin, zmax, npart_des, npart_out, stretch, thres, pos, pvol, nu, h, get_density, get_id, validate_position)

    Places particles on a lattice containing a physical object

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    double precision, intent(in) :: central_density

    Maximum baryon mass density of the system

    double precision, intent(in) :: xmin

    Left boundary of the lattice

    double precision, intent(in) :: xmax

    Right boundary of the lattice

    double precision, intent(in) :: ymin

    Left boundary of the lattice

    double precision, intent(in) :: ymax

    Right boundary of the lattice

    double precision, intent(in) :: zmin

    Left boundary of the lattice

    double precision, intent(in) :: zmax

    Right boundary of the lattice

    integer, intent(in) :: npart_des

    Desired particle number

    integer, intent(out) :: npart_out

    Real, output particle number

    double precision, intent(in) :: stretch

    Stretching factor fo the lattice. xmin to zmax are multiplied by it

    double precision, intent(in) :: thres

    (~rho_max)/thres is the minimum mass density considered when placing particles. Used only when redistribute_nu is .FALSE. . When redistribute_nu is .TRUE. thres= 100*nu_ratio

    double precision, intent(inout), DIMENSION(:,:), ALLOCATABLE :: pos
    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: pvol

    Array storing the particle positions Array storing the particle volumes

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: nu

    Array storing the particle baryon masses

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: h

    Array storing the initial guess for the particle smoothing lengths

    function get_density(x, y, z) result(density)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    Return Value double precision

    Baryon mass density at

    subroutine get_id(x, y, z, sqdetg, baryon_density, gamma_euler)

    Returns the baryon mass density at the desired point

    Arguments
    Type IntentOptional Attributes Name
    double precision, intent(in) :: x

    coordinate of the desired point

    double precision, intent(in) :: y

    coordinate of the desired point

    double precision, intent(in) :: z

    coordinate of the desired point

    double precision, intent(out) :: sqdetg
    double precision, intent(out) :: baryon_density
    double precision, intent(out) :: gamma_euler
    procedure(validate_position_int), optional :: validate_position

    Returns 1 if the position is not valid, 0 otherwise

interface

  • public module subroutine print_formatted_id_particles(this, namefile)

    Prints the SPH ID to a formatted file

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    character(len=*), intent(inout), optional :: namefile

    Name of the formatted output file

interface

  • public module subroutine print_summary(this)

    Prints a summary of the properties of the \(\mathrm{SPH}\) particle distribution, optionally, to a formatted file whose name is given as the optional argument filename

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    Name of the formatted file to print the summary to

interface

  • public module subroutine read_compose_composition(this, namefile)

    Reads the table in the CompOSE file with extension .beta

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    character(len=*), intent(inout), optional :: namefile

    To read the file great_eos.beta in directory compose_path/GREAT_EoS, namefile="GREAT_EoS/great_eos"

interface

  • public module subroutine read_particles_formatted_file(this, parts_pos_unit, nline_in, nline_fin, xmin, xmax, ymin, ymax, zmin, zmax, pos, pvol, nu, h)

    Read particle positions and nu from a formatted file with the following format:

    1. The first line must contain the integer number of objects for the system, and the particle numbers on each object; each column separated by one tab.
    2. The other lines must contain at least 3 columns with the x, y, z coordinates of the particle positions. An optional 4th column can contain the values of nu on the particles. If the 4th column is not present, nu will be roughly estimated using the density read from the ID and the average volume per particle; the latter is computed by computing the average particle separation along the z direction.

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    String containing the name of the formatted file containing the particle positions and, optionally, nu

    integer :: parts_pos_unit

    Unit number of the formatted file containing the particle positions and, optionally, nu

    integer :: nline_in

    First line containing the relevant data

    integer :: nline_fin

    Last line containing the relevant data

    double precision, intent(in) :: xmin

    Left boundary of the lattice

    double precision, intent(in) :: xmax

    Right boundary of the lattice

    double precision, intent(in) :: ymin

    Left boundary of the lattice

    double precision, intent(in) :: ymax

    Right boundary of the lattice

    double precision, intent(in) :: zmin

    Left boundary of the lattice

    double precision, intent(in) :: zmax

    Right boundary of the lattice

    double precision, intent(inout), DIMENSION(:,:), ALLOCATABLE :: pos
    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: pvol

    Array storing the particle positions Array storing the particle volumes

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: nu

    Array storing the particle baryon masses

    double precision, intent(inout), DIMENSION(:), ALLOCATABLE :: h

    Array storing the initial guess for the particle smoothing lengths

interface

  • public module subroutine read_sphincs_dump_print_formatted(this, namefile_bin, namefile, save_data)

    Reads the binary ID file printed by compute_and_print_sph_variables and prints the data stored in it to a formatted file

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    character(len=*), intent(inout), optional :: namefile_bin

    Name of the binary file to be read

    character(len=*), intent(inout), optional :: namefile

    Name of the formatted file to be printed

    logical, intent(in), optional :: save_data

    If .TRUE., saves the read data into the TYPE member variables

interface

  • public module subroutine reflect_particles_xy_plane(npart, pos, pos_below, npart_above_xy, above_xy_plane_a, nu, nu_below)

    Mirror the particle with z>0 with respect to the xy plane, to impose the equatorial-plane symmetry

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart
    double precision, intent(inout), DIMENSION(3,npart) :: pos
    double precision, intent(out), DIMENSION(3,npart_above_xy) :: pos_below
    integer, intent(inout) :: npart_above_xy
    integer, intent(in), DIMENSION(npart_above_xy) :: above_xy_plane_a
    double precision, intent(inout), optional, DIMENSION(npart) :: nu
    double precision, intent(out), optional, DIMENSION(npart_above_xy) :: nu_below

interface

  • public module subroutine test_recovery(this, npart, pos, nlrf, u, pr, vel_u, theta, nstar)

    Tests the recovery. Computes the conserved variables from the physical ones, and then the physical ones from the conserved ones. It then compares the variables computed with the recovery PROCEDURES, with those computed with \(\texttt{SPHINCS_ID}\).

    Arguments

    Type IntentOptional Attributes Name
    class(particles), intent(inout) :: this

    particles object which this PROCEDURE is a member of

    integer, intent(in) :: npart

    Particle number

    double precision, intent(in), DIMENSION(3,npart) :: pos

    Particle positions

    double precision, intent(in), DIMENSION(npart) :: nlrf

    Baryon density in the local rest frame on the particles

    double precision, intent(in), DIMENSION(npart) :: u

    Specific internal energy on the particles

    double precision, intent(in), DIMENSION(npart) :: pr

    Pressure on the particles

    double precision, intent(in), DIMENSION(3,npart) :: vel_u

    Spatial velocity in the computing frame on the particles

    double precision, intent(in), DIMENSION(npart) :: theta

    Generalized Lorentz factor on the particles

    double precision, intent(in), DIMENSION(npart) :: nstar

    Proper baryon density in the local rest frame on the particles

    Canonical momentum on the particles

    Canonical energy on the particles

    Name of the formatted file where the data is printed


Abstract Interfaces

abstract interface

  • public subroutine post_process_sph_id_int(npart, pos, nlrf, u, pr, vel_u, theta, nstar, nu)

    Post-process the \(\mathrm{SPH}\) \(\mathrm{ID}\); for example, correct for the residual ADM linear momentum.

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: npart

    Particle number

    double precision, intent(inout), DIMENSION(3,npart) :: pos

    Particle positions

    double precision, intent(inout), DIMENSION(npart) :: nlrf

    Baryon density in the local rest frame on the particles

    double precision, intent(inout), DIMENSION(npart) :: u

    Specific internal energy on the particles

    double precision, intent(inout), DIMENSION(npart) :: pr

    Pressure on the particles

    double precision, intent(inout), DIMENSION(3,npart) :: vel_u

    Spatial velocity in the computing frame on the particles

    double precision, intent(inout), DIMENSION(npart) :: theta

    Generalized Lorentz factor on the particles

    double precision, intent(inout), DIMENSION(npart) :: nstar

    Proper baryon density in the local rest frame on the particles

    double precision, intent(inout), DIMENSION(npart) :: nu

    Baryon number per particle


Derived Types

type, public ::  eos

Data structure representing an \(\mathrm{EOS}\)

Components

Type Visibility Attributes Name Initial
character(len=:), public, ALLOCATABLE :: eos_name

The \(\mathrm{EOS}\) name

double precision, public, DIMENSION(:), ALLOCATABLE :: eos_parameters

Array containing the \(\mathrm{EOS}\) parameters, in the following order:

Read more…
double precision, public, DIMENSION(:,:), ALLOCATABLE :: table_eos

Array containing the tabulated \(\mathrm{EOS}\). Its size is determined when reading the tabulated \(\mathrm{EOS}\) from file. The first index runs over the field, the second over its values.

type, public ::  particles

TYPE representing a \(\mathrm{SPH}\) particle distribution

Components

Type Visibility Attributes Name Initial
double precision, private, DIMENSION(:), ALLOCATABLE :: Theta

1-D array storing the generalized Lorentz factor

double precision, private, DIMENSION(:), ALLOCATABLE :: Ye

1-D array storing the electron fraction

double precision, private, DIMENSION(:), ALLOCATABLE :: Ye_table

Array storing the values of the electron fraction in the \(\mathrm{CompOSE}\) table

double precision, private, DIMENSION(3) :: adm_linear_momentum_fluid

Estimate of the linear momentum of the fluid computed from the canonical momentum per baryon on the particles

double precision, private, DIMENSION(:,:), ALLOCATABLE :: adm_linear_momentum_i

Estimate of the linear momentum of each matter object, computed from the canonical momentum per baryon on the particles

double precision, private :: adm_mass
type(eos), private, DIMENSION(:), ALLOCATABLE :: all_eos

Array of TYPE eos containing the \(\mathrm{EOS}\) information for all the matter objects

logical, private, DIMENSION(:), ALLOCATABLE :: apm_iterate

.TRUE. if the Artificial Pressure Method (APM) has to be applied to the particles, .FALSE. otherwise

type(timer), public, DIMENSION(:), ALLOCATABLE :: apm_timers

Timer that times how long it takes to perform the APM on the matter objects

double precision, private, DIMENSION(:,:), ALLOCATABLE :: barycenter

Array storing the centers of mass of the matter objects

double precision, private, DIMENSION(4) :: barycenter_system

Array storing the center of mass of the entire particle distribution

double precision, private, DIMENSION(:), ALLOCATABLE :: baryon_density

1-D array storing the baryon mass density in the fluid frame

integer, private, DIMENSION(:), ALLOCATABLE :: baryon_density_index

Array storing the indices to use with baryon_density to sort the elements of baryon_density in increasing order

integer, private :: call_flag = 0

Flag that is set different than 0 if the SUBROUTINE compute_and_print_sph_variables is called

logical, private :: cold_system

.TRUE. if the system is at zero temperature (no thermal component); .FALSE. otherwise

logical, private :: compose_eos

.TRUE. if the electron fraction should be read from the CompOSE table with extension.beta, .FALSE. otherwise

Read more…
character(len=:), private, ALLOCATABLE :: compose_filename

String storing the subpath of compose_path to the CompOSE file with .beta extension

character(len=:), private, ALLOCATABLE :: compose_path

String storing the local path to the directory containing the CompOSE \(\mathrm{EOS}\)

logical, private :: correct_nu

.TRUE. if the baryon number per particle should be corrected to account for the total baryon masses of the stars, .FALSE. otherwise

integer, private :: distribution_id

Identification number for the particle distribution

logical, private :: empty_object

.TRUE. if the object is empty, .FALSE. if it's not empty

double precision, private, DIMENSION(:), ALLOCATABLE :: energy_density

1-D array storing the energy density

double precision, private, DIMENSION(:), ALLOCATABLE :: enthalpy

1-D array storing the enthalpy

Read more…
logical, public :: export_bin

.TRUE. if the binary files for SPHINCS_BSSN are to be exported, .FALSE. otherwise

logical, public :: export_form_x

.TRUE. if the ID in the formatted files is to be on the x axis only, .FALSE. otherwise

logical, public :: export_form_xy

.TRUE. if the ID in the formatted files is to be on the xy plane only, .FALSE. otherwise

double precision, private, DIMENSION(:), ALLOCATABLE :: g_xx

Array storing the values of the xx component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: g_xy

Array storing the values of the xy component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: g_xz

Array storing the values of the xz component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: g_yy

Array storing the values of the xz component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: g_yz

Array storing the values of the yz component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: g_zz

Array storing the values of the zz component of the spatial metric on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: h

1-D array storing the smoothing length

type(timer), public :: importer_timer

Timer that times how long it takes to import the ID at the particle positions

double precision, private, DIMENSION(:), ALLOCATABLE :: lapse

Array storing the values of the lapse function on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: mass_fractions

Ratio between baryonic masses of the matter objects and the total baryonic mass of all the matter objects

Read more…
double precision, private, DIMENSION(:), ALLOCATABLE :: mass_ratios

Ratio between baryonic masses of the matter objects and the maximum baryonic mass among them

Read more…
double precision, private, DIMENSION(:), ALLOCATABLE :: masses

1-D array storing the particle masses

integer, private :: n_matter

Number of matter objects in the physical system

double precision, private, DIMENSION(:), ALLOCATABLE :: nb_table

Array storing the values of the baryon number density in the \(\mathrm{CompOSE}\) table.

Read more…
double precision, private, DIMENSION(:), ALLOCATABLE :: nbar_i

Baryon numbers of the matter objects

double precision, private :: nbar_tot

Total baryon number

double precision, private, DIMENSION(:), ALLOCATABLE :: nlrf

1-D array storing baryon density in the local rest frame , computed directly from the \(\texttt{LORENE}\) density

double precision, private, DIMENSION(:), ALLOCATABLE :: nlrf_sph

1-D array storing baryon density in the local rest frame , computed from the kernel interpolated proper baryon number density nstar_sph

integer, private :: npart

Total particle number

integer, private, DIMENSION(:), ALLOCATABLE :: npart_fin

Array storing the index of the last particle for each matter objects

integer, private, DIMENSION(:), ALLOCATABLE :: npart_i

Array storing the particle numbers for the matter objects

double precision, private, DIMENSION(:), ALLOCATABLE :: nstar

1-D array storing the SPH estimate of the proper baryon number density

double precision, private, DIMENSION(:), ALLOCATABLE :: nstar_sph

1-D array storing the SPH estimate of the proper baryon number density, from kernel interpolation

double precision, private, DIMENSION(:), ALLOCATABLE :: nu

1-D array storing the baryon number per particle

double precision, private :: nu_ratio_des

Desired ratio between the max and min of the baryon number per particle, over all the matter objects. Only used when redistribute_nu is .TRUE.

Read more…
double precision, private :: nuratio

Ratio between the max and min of the baryon number per particle, over all the matter objects

double precision, private, DIMENSION(:), ALLOCATABLE :: nuratio_i

Baryon number ratios on the matter objects

double precision, private, DIMENSION(:), ALLOCATABLE :: particle_density

1-D array storing the particle number density

double precision, private, DIMENSION(:), ALLOCATABLE :: particle_density_sph

1-D array storing the SPH estimate of the particle number density, from kernel interpolation

type(timer), public :: placer_timer

Timer that times how long it takes to place particles on the stars

double precision, private, DIMENSION(:,:), ALLOCATABLE :: pos

2-D array storing the particle positions

procedure, private, POINTER, NOPASS :: post_process_sph_id

Pointer to a procedure that post_process the \(\mathrm{SPH}\) \(\mathrm{ID}\); for example, correct for the residual ADM linear momentum.

double precision, private, DIMENSION(:), ALLOCATABLE :: pressure

1-D array storing the pressure

double precision, private, DIMENSION(:), ALLOCATABLE :: pressure_sph

1-D array storing the pressure in code units

double precision, private, DIMENSION(:), ALLOCATABLE :: pvol

1-D array storing the particle volumes

logical, private :: randomize_phi

.TRUE. if the particle positions on ellipsoidal surfaces are randomized in the direction, .FALSE. otherwise

logical, private :: randomize_r

.TRUE. if the particle positions on ellipsoidal surfaces are randomized in the direction, .FALSE. otherwise

logical, private :: randomize_theta

.TRUE. if the particle positions on ellipsoidal surfaces are randomized in the direction, .FALSE. otherwise

logical, private :: read_nu

.TRUE. if the baryon number per particle has to be read from the formatted file containing the particle positions, .FALSE. otherwise

logical, private :: redistribute_nu

.TRUE. if the baryon number per particle should be reassigned, trying to obtain a baryon number ratio no larger than nu_ratio, when placing particles on lattices; .FALSE. otherwise

logical, private :: reflect_particles_x

.TRUE. if the particles on star 2 should be the reflection of the particles on star 1 with respect to the plane, only if the baryon masses of the stars differe less than ; .FALSE. otherwise

type(timer), public :: same_particle_timer

Timer that times how long it takes to check if there are multiple particles at the same positions

double precision, private, DIMENSION(:), ALLOCATABLE :: shift_x

Array storing the values of the x component of the shift vector on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: shift_y

Array storing the values of the y component of the shift vector on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: shift_z

Array storing the values of the z component of the shift vector on the particles

double precision, private, DIMENSION(:), ALLOCATABLE :: specific_energy

1-D array storing the specific internal energy

type(timer), public :: sph_computer_timer

Timer that times how long it takes to compute the SPH variables at the particle positions

character(len=50), private :: sphincs_id_particles

String containing the name of the particles parameter file

type(surface), private, DIMENSION(:), ALLOCATABLE :: surfaces

1-D array storing the surfaces of the matter objects

double precision, private, DIMENSION(:), ALLOCATABLE :: u_sph

1-D array storing the specific internal energy computed using formula (9) in Read et al., Phys.Rev.D79:124032,2009, [arXiv:0812.2163][https://arxiv.org/abs/0812.2163]{:target="_blank"}

logical, private, DIMENSION(:), ALLOCATABLE :: use_atmosphere

.TRUE. to allow the particles to move where the density is 0 during the Artificial Pressure Method (APM) iteration. This can be useful when the system has an irregular geometry, as, for example, an ejecta .FALSE. otherwise

logical, private :: use_thres

.TRUE. if the threshold on the baryon mass density should e applied when placing particles on lattices, .FALSE. otherwise

double precision, private, DIMENSION(:,:), ALLOCATABLE :: v

2-D array storing the coordinate fluid 4-velocity

double precision, private, DIMENSION(:), ALLOCATABLE :: v_euler_x

1-D array storing the x component of the fluid 3-velocity wrt the Eulerian observer

double precision, private, DIMENSION(:), ALLOCATABLE :: v_euler_y

1-D array storing the y component of the fluid 3-velocity wrt the Eulerian observer

double precision, private, DIMENSION(:), ALLOCATABLE :: v_euler_z

1-D array storing the z component of the fluid 3-velocity wrt the Eulerian observer

Constructor

Interface of TYPE particles

public interface construct_particles_std ()
public interface construct_particles_bin ()

Finalizations Procedures

final :: destruct_particles

Finalizer (Destructor) of particles object

Type-Bound Procedures

procedure , public :: allocate_particles_memory Interface

Allocates memory for the particles member arrays

procedure , public :: analyze_hydro Interface

Scans the hydro fields taken from to look for negative or zero values

procedure , public :: compute_Ye Interface

Interpates linearly the electron fraction at the particle densities; that is, assigns at the particle positions

procedure , public :: compute_and_print_sph_variables Interface

Computes the SPH variables at the particle positions, and optionally prints them to a binary file to be read by and , and to a formatted file to be read by , by calling print_formatted_id_particles

procedure , public :: compute_sph_hydro Interface

Computes the hydro fields on a section of the particles specified as input. First, computes the \(\mathrm{SPH}\) pressure starting from the \(\mathrm{SPH}\) baryon mass density, and the specific internal energy. The pressure is computed differently for different \(\mathrm{EOS}\), and for cold and hot systems. Then computes the enthalpy and the sound speed accordingly.

procedure , public :: deallocate_particles_memory Interface

Deallocates memory for the particles member arrays

procedure , public :: get_compose_eos Interface

Returns compose_eos

procedure , public :: get_eos_id Interface

Returns the \(\mathrm{EOS}\) identifier for matter object i_matter

procedure , public :: get_g3 Interface

Returns (g_xx,g_xy,g_xz,g_yy,g_yz,g_zz)

procedure , public :: get_h Interface

Returns h

procedure , public :: get_lapse Interface

Returns lapse

procedure , public :: get_n_matter Interface

Returns n_matter

procedure , public :: get_nlrf Interface

Returns nlrf

procedure , public :: get_nlrf_sph Interface

Returns nlrf_sph

procedure , public :: get_npart Interface

Returns npart

procedure , public :: get_npart_i Interface

Returns the number of particles on the object i_matter

procedure , public :: get_nstar Interface

Returns nstar

procedure , public :: get_nstar_sph Interface

Returns nstar_sph

procedure , public :: get_nu Interface

Returns nu

procedure , public :: get_nuratio Interface

Returns nuratio

procedure , public :: get_nuratio_i Interface

Returns the baryon number ratio on the object i_matter

procedure , public :: get_object_of_particle Interface

Returns the number of the matter object asociated with the particle number given as input. Example: give number as input; this particle number corresponds to a particle on matter object . This functions returns .

procedure , public :: get_pos Interface

Returns pos

procedure , public :: get_pressure Interface

Returns pressure

procedure , public :: get_pressure_sph Interface

Returns pressure_sph

procedure , public :: get_shift Interface

Returns (shift_x,shift_y,shift_z)

procedure , public :: get_theta Interface

Returns Theta

procedure , public :: get_u Interface

Returns specific_energy

procedure , public :: get_u_sph Interface

Returns u_sph

procedure , public :: get_vel Interface

Returns v

procedure , public , NOPASS :: impose_equatorial_plane_symmetry Interface

Reflects the positions of the particles on a matter object with respect to the plane

procedure , public :: is_empty Interface

Returns .TRUE if the particles object is empty, .FALSE otherwise

Read more…
procedure , public , NOPASS :: perform_apm Interface

Performs the Artificial Pressure Method (APM) on one star's particles

procedure , public :: place_particles_ellipsoidal_surfaces Interface

Places particles on ellipsoidal surfaces on one star

procedure , public :: place_particles_lattice Interface

Places particles on a single lattice that surrounds both stars

procedure , public :: print_formatted_id_particles Interface

Prints the \(\mathrm{SPH}\) \(\mathrm{ID}\) to a formatted file

procedure , public :: print_summary Interface

Prints the SPH ID to a formatted file

procedure , public :: read_compose_composition Interface

Reads the table in the CompOSE file with extension .beta

procedure , public :: read_particles_formatted_file Interface

Read particle positions and nu from a formatted file with the following format:

Read more…
procedure , public :: read_sphincs_dump_print_formatted Interface

Reads the binary ID file printed by compute_and_print_sph_variables, and prints it to a formatted file.

Read more…
procedure , public :: test_recovery Interface

Computes the conserved variables from the physical ones, and vice versa, to test that the recovered physical variables are the same to those computed from the \(\mathrm{ID}\).

Read more…

Functions

public function find_h_backup(a, npart, pos, ndes) result(h)

Backup method to find the smoothing length via brute force if the optimized method gives 0. It sets the smoothing lengths to the distance between the particle and the ndes-th closest neighbour.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: a

Index of the particle whose smoothing length is to be computed

integer, intent(in) :: npart

Number of particles

double precision, intent(in), DIMENSION(3,npart) :: pos

Array containing particle positions

integer, intent(in) :: ndes

Desired number of neighbours

Return Value double precision

Smoothing length