construct_particles_std Module Procedure

module procedure construct_particles_std module function construct_particles_std(id, dist) result(parts)

Uses

    • options
    • units
    • utility
    • input_output
    • constants
    • alive_flag
    • analyze
    • kernel_table
  • proc~~construct_particles_std~~UsesGraph proc~construct_particles_std construct_particles_std alive_flag alive_flag proc~construct_particles_std->alive_flag analyze analyze proc~construct_particles_std->analyze constants constants proc~construct_particles_std->constants input_output input_output proc~construct_particles_std->input_output kernel_table kernel_table proc~construct_particles_std->kernel_table module~utility utility proc~construct_particles_std->module~utility options options proc~construct_particles_std->options units units proc~construct_particles_std->units module~utility->constants matrix matrix module~utility->matrix

The constructor of TYPE particles is supposed to set up a particle distribution by assigning the particle positions, their baryon numbers nu and first guesses for their smoothing lengths h. It also sets up the unit system and the kernel.

After the particle distribution is set up, it assigns the \(\mathrm{ID}\) to the particles. It does not compute the \(\mathrm{SPH}\) variables and it does not set up the neighbors' tree. The latter two things are delegated to the other methods of TYPE particles.

FT 17.10.2020


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


Calls

proc~~construct_particles_std~~CallsGraph proc~construct_particles_std construct_particles_std alive alive proc~construct_particles_std->alive com com proc~construct_particles_std->com interface~check_particle_positions check_particle_positions proc~construct_particles_std->interface~check_particle_positions interface~perform_apm perform_apm proc~construct_particles_std->interface~perform_apm none~check_eos check_eos proc~construct_particles_std->none~check_eos none~place_particles_on_ellipsoidal_surfaces place_particles_on_ellipsoidal_surfaces proc~construct_particles_std->none~place_particles_on_ellipsoidal_surfaces none~place_particles_on_lattices place_particles_on_lattices proc~construct_particles_std->none~place_particles_on_lattices none~read_particles_from_formatted_file read_particles_from_formatted_file proc~construct_particles_std->none~read_particles_from_formatted_file none~read_particles_options read_particles_options proc~construct_particles_std->none~read_particles_options none~reflect_particles_yz_plane reflect_particles_yz_plane proc~construct_particles_std->none~reflect_particles_yz_plane proc~scan_1d_array_for_nans scan_1d_array_for_nans proc~construct_particles_std->proc~scan_1d_array_for_nans proc~spatial_vector_norm_sym3x3 spatial_vector_norm_sym3x3 proc~construct_particles_std->proc~spatial_vector_norm_sym3x3 proc~spherical_from_cartesian spherical_from_cartesian proc~construct_particles_std->proc~spherical_from_cartesian timer timer proc~construct_particles_std->timer proc~check_particle_positions check_particle_positions interface~check_particle_positions->proc~check_particle_positions proc~perform_apm perform_apm interface~perform_apm->proc~perform_apm none~place_particles_on_ellipsoidal_surfaces->none~reflect_particles_yz_plane none~place_particles_on_lattices->none~reflect_particles_yz_plane none~read_particles_from_formatted_file->none~reflect_particles_yz_plane interface~impose_equatorial_plane_symmetry impose_equatorial_plane_symmetry none~read_particles_from_formatted_file->interface~impose_equatorial_plane_symmetry ktable ktable none~read_particles_options->ktable read_options read_options none~read_particles_options->read_options set_units set_units none~read_particles_options->set_units proc~is_finite_number is_finite_number proc~scan_1d_array_for_nans->proc~is_finite_number proc~impose_equatorial_plane_symmetry impose_equatorial_plane_symmetry interface~impose_equatorial_plane_symmetry->proc~impose_equatorial_plane_symmetry indexx indexx proc~check_particle_positions->indexx proc~perform_apm->com proc~perform_apm->interface~check_particle_positions proc~perform_apm->proc~spherical_from_cartesian proc~perform_apm->timer proc~perform_apm->interface~impose_equatorial_plane_symmetry proc~perform_apm->proc~is_finite_number allocate_gradient allocate_gradient proc~perform_apm->allocate_gradient allocate_metric_on_particles allocate_metric_on_particles proc~perform_apm->allocate_metric_on_particles allocate_rcb_tree_memory_3d allocate_rcb_tree_memory_3d proc~perform_apm->allocate_rcb_tree_memory_3d allocate_sph_memory allocate_sph_memory proc~perform_apm->allocate_sph_memory assign_h assign_h proc~perform_apm->assign_h deallocate_gradient deallocate_gradient proc~perform_apm->deallocate_gradient deallocate_metric_on_particles deallocate_metric_on_particles proc~perform_apm->deallocate_metric_on_particles deallocate_rcb_tree_memory_3d deallocate_rcb_tree_memory_3d proc~perform_apm->deallocate_rcb_tree_memory_3d deallocate_sph_memory deallocate_sph_memory proc~perform_apm->deallocate_sph_memory density density proc~perform_apm->density density_loop density_loop proc~perform_apm->density_loop exact_nei_tree_update exact_nei_tree_update proc~perform_apm->exact_nei_tree_update h h proc~perform_apm->h interface~correct_center_of_mass correct_center_of_mass proc~perform_apm->interface~correct_center_of_mass iorig iorig proc~perform_apm->iorig none~allocate_apm_fields allocate_apm_fields proc~perform_apm->none~allocate_apm_fields none~discard_atmosphere discard_atmosphere proc~perform_apm->none~discard_atmosphere none~dump_apm_pos dump_apm_pos proc~perform_apm->none~dump_apm_pos none~get_nstar_id_atm get_nstar_id_atm proc~perform_apm->none~get_nstar_id_atm none~place_and_print_ghost_particles place_and_print_ghost_particles proc~perform_apm->none~place_and_print_ghost_particles none~read_pressure_id read_pressure_id proc~perform_apm->none~read_pressure_id none~reallocate_output_fields reallocate_output_fields proc~perform_apm->none~reallocate_output_fields none~validate_position_final validate_position_final proc~perform_apm->none~validate_position_final nu nu proc~perform_apm->nu position_correction position_correction proc~perform_apm->position_correction proc~cartesian_from_spherical cartesian_from_spherical proc~perform_apm->proc~cartesian_from_spherical proc~find_h_backup find_h_backup proc~perform_apm->proc~find_h_backup proc~correct_center_of_mass correct_center_of_mass interface~correct_center_of_mass->proc~correct_center_of_mass none~discard_atmosphere->h none~discard_atmosphere->nu get_density get_density none~discard_atmosphere->get_density none~dump_apm_pos->get_density none~get_nstar_id_atm->proc~is_finite_number get_nstar_id get_nstar_id none~get_nstar_id_atm->get_nstar_id none~place_and_print_ghost_particles->proc~spherical_from_cartesian none~place_and_print_ghost_particles->timer none~place_and_print_ghost_particles->proc~is_finite_number none~place_and_print_ghost_particles->allocate_gradient none~place_and_print_ghost_particles->allocate_metric_on_particles none~place_and_print_ghost_particles->allocate_rcb_tree_memory_3d none~place_and_print_ghost_particles->allocate_sph_memory none~place_and_print_ghost_particles->assign_h none~place_and_print_ghost_particles->deallocate_gradient none~place_and_print_ghost_particles->deallocate_metric_on_particles none~place_and_print_ghost_particles->deallocate_rcb_tree_memory_3d none~place_and_print_ghost_particles->deallocate_sph_memory none~place_and_print_ghost_particles->density_loop none~place_and_print_ghost_particles->h none~place_and_print_ghost_particles->iorig none~place_and_print_ghost_particles->proc~cartesian_from_spherical none~place_and_print_ghost_particles->proc~find_h_backup bilinear_interpolation bilinear_interpolation none~place_and_print_ghost_particles->bilinear_interpolation center center none~place_and_print_ghost_particles->center compute_pressure compute_pressure none~place_and_print_ghost_particles->compute_pressure none~place_and_print_ghost_particles->get_density none~place_and_print_ghost_particles->get_nstar_id nu_output nu_output none~place_and_print_ghost_particles->nu_output pos_input pos_input none~place_and_print_ghost_particles->pos_input sizes sizes none~place_and_print_ghost_particles->sizes get_pressure_id get_pressure_id none~read_pressure_id->get_pressure_id h_output h_output none~reallocate_output_fields->h_output none~reallocate_output_fields->nu_output none~reallocate_output_fields->pos_input validate_position validate_position none~validate_position_final->validate_position proc~impose_equatorial_plane_symmetry->com interface~find_particles_above_xy_plane find_particles_above_xy_plane proc~impose_equatorial_plane_symmetry->interface~find_particles_above_xy_plane interface~reflect_particles_xy_plane reflect_particles_xy_plane proc~impose_equatorial_plane_symmetry->interface~reflect_particles_xy_plane proc~find_particles_above_xy_plane find_particles_above_xy_plane interface~find_particles_above_xy_plane->proc~find_particles_above_xy_plane proc~reflect_particles_xy_plane reflect_particles_xy_plane interface~reflect_particles_xy_plane->proc~reflect_particles_xy_plane proc~correct_center_of_mass->com proc~correct_center_of_mass->proc~is_finite_number

Called by

proc~~construct_particles_std~~CalledByGraph proc~construct_particles_std construct_particles_std interface~construct_particles_std construct_particles_std interface~construct_particles_std->proc~construct_particles_std interface~particles particles interface~particles->interface~construct_particles_std program~convergence_test convergence_test program~convergence_test->interface~particles program~sphincs_id sphincs_id program~sphincs_id->interface~particles

Contents


Variables

Type Visibility Attributes Name Initial
integer, private :: a
logical, private :: adapt_ghosts
logical, private, DIMENSION(max_length) :: apm_iterate
integer, private :: apm_max_it
double precision, private, DIMENSION(id% get_n_matter(),3) :: barycenter
double precision, private, DIMENSION(id% get_n_matter(),3) :: center
double precision, private, DIMENSION(id% get_n_matter()) :: central_density
integer, private :: column_nu
integer, private, DIMENSION(3) :: columns
logical, private :: compose_eos
character(len=max_length), private :: compose_filename
character(len=max_length), private :: compose_path
logical, private :: correct_nu
integer, private, SAVE :: counter = 1
logical, private, parameter :: debug = .FALSE.
logical, private :: exist
logical, private :: file_exists
character(len=max_length), private :: filename_apm_pos
character(len=max_length), private :: filename_apm_pos_id
character(len=max_length), private :: filename_apm_results
character(len=max_length), private :: filename_mass_profile
character(len=max_length), private :: filename_shells_pos
character(len=max_length), private :: filename_shells_radii
double precision, private, DIMENSION(id% get_n_matter()) :: g00_lengthscales
double precision, private :: ghost_dist
double precision, private, DIMENSION(id% get_n_matter()) :: ghost_dists
integer, private :: header_lines
integer, private :: i_matter
double precision, private, DIMENSION(id% get_n_matter()) :: lapse_lengthscales
double precision, private :: last_r
double precision, private :: lower_bound
double precision, private :: lower_factor
logical, private :: mass_it
integer, private :: max_inc
integer, private, parameter :: max_length = 50
double precision, private :: max_mass
integer, private :: max_steps
double precision, private :: min_eps
double precision, private :: min_g00_abs
double precision, private :: min_lapse
double precision, private :: min_vel
logical, private :: move_away_ghosts
integer, private :: n_cols
integer, private :: n_matter_tmp
integer, private :: nlines
integer, private :: npart_des
integer, private, DIMENSION(id% get_n_matter()) :: npart_des_i
integer, private, DIMENSION(:), ALLOCATABLE :: npart_i_tmp
integer, private :: npart_tmp
double precision, private :: nu_ratio_des
double precision, private :: nuratio_des
double precision, private :: nuratio_thres
integer, private :: nx_gh
integer, private :: ny_gh
integer, private :: nz_gh
type(parts_i), private, DIMENSION(id% get_n_matter()) :: parts_all
character(len=:), private, ALLOCATABLE :: parts_out_namefile

Name for the file to print the final particle distribution and nu

character(len=max_length), private :: parts_pos
character(len=:), private, ALLOCATABLE :: parts_pos_namefile
character(len=max_length), private :: parts_pos_path
double precision, private :: phi_a
integer, private :: print_step
double precision, private :: r_a
logical, private :: randomize_phi
logical, private :: randomize_r
logical, private :: randomize_theta
logical, private :: read_nu
logical, private :: redistribute_nu
logical, private :: reflect_particles_x
logical, private, DIMENSION(max_length) :: remove_atmosphere
double precision, private :: shift_norm
double precision, private, DIMENSION(id% get_n_matter(),6) :: sizes
character(len=3), private :: str_i
double precision, private :: stretch
double precision, private :: theta_a
double precision, private :: thres
integer, private :: tmp
double precision, private, parameter :: tol_equal_mass = 5.0D-3
double precision, private :: total_mass
integer, private, parameter :: unit_pos = 2289
integer, private, parameter :: unit_pos_out = 8754
double precision, private :: upper_bound
double precision, private :: upper_factor
logical, private, DIMENSION(max_length) :: use_atmosphere
logical, private :: use_pressure
logical, private :: use_thres
double precision, private :: xmax
double precision, private :: xmin
double precision, private :: ymax
double precision, private :: ymin
double precision, private :: zmax
double precision, private :: zmin

Functions

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

Wrapper function to read the baryon mass density from the \(\mathrm{ID}\)

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 import_pressure_id(x, y, z) result(pressure)

Wrapper function to read the pressure from the \(\mathrm{ID}\)

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_position(x, y, z) result(answer)

Wrapper function to validate a position

Arguments

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

Return Value logical


Subroutines

subroutine check_eos()

Check that the supplied \(\mathrm{EOS}\) parameters are consistent with the \(\mathrm{EOS}\) used to compute the \(\mathrm{ID}\)

FT xx.09.2022


Arguments

None

subroutine compute_nstar_eul_id(npart, v_euler_x, v_euler_y, v_euler_z, g_xx, g_xy, g_xz, g_yy, g_yz, g_zz, baryon_density, nstar_eul_id)

Compute nstar_eul_id, the relativistic baryon mass density seen by the Eulerian observer, given the \(\mathrm{ID}\)

FT 31.08.2021


Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npart
double precision, intent(in), DIMENSION(npart) :: v_euler_x
double precision, intent(in), DIMENSION(npart) :: v_euler_y
double precision, intent(in), DIMENSION(npart) :: v_euler_z
double precision, intent(in), DIMENSION(npart) :: g_xx
double precision, intent(in), DIMENSION(npart) :: g_xy
double precision, intent(in), DIMENSION(npart) :: g_xz
double precision, intent(in), DIMENSION(npart) :: g_yy
double precision, intent(in), DIMENSION(npart) :: g_yz
double precision, intent(in), DIMENSION(npart) :: g_zz
double precision, intent(in), DIMENSION(npart) :: baryon_density
double precision, intent(out), DIMENSION(npart) :: nstar_eul_id

subroutine compute_nstar_id(npart, lapse, shift_x, shift_y, shift_z, v_euler_x, v_euler_y, v_euler_z, g_xx, g_xy, g_xz, g_yy, g_yz, g_zz, baryon_density, nstar_sph, nstar_id, nlrf_sph, sqg)

Compute nstar_id, the relativistic baryon mass density, given the required \(\mathrm{ID}\) as input

FT 31.08.2021


Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npart
double precision, intent(in), DIMENSION(npart) :: lapse
double precision, intent(in), DIMENSION(npart) :: shift_x
double precision, intent(in), DIMENSION(npart) :: shift_y
double precision, intent(in), DIMENSION(npart) :: shift_z
double precision, intent(in), DIMENSION(npart) :: v_euler_x
double precision, intent(in), DIMENSION(npart) :: v_euler_y
double precision, intent(in), DIMENSION(npart) :: v_euler_z
double precision, intent(in), DIMENSION(npart) :: g_xx
double precision, intent(in), DIMENSION(npart) :: g_xy
double precision, intent(in), DIMENSION(npart) :: g_xz
double precision, intent(in), DIMENSION(npart) :: g_yy
double precision, intent(in), DIMENSION(npart) :: g_yz
double precision, intent(in), DIMENSION(npart) :: g_zz
double precision, intent(in), DIMENSION(npart) :: baryon_density
double precision, intent(in), DIMENSION(npart) :: nstar_sph
double precision, intent(out), DIMENSION(npart) :: nstar_id
double precision, intent(out), DIMENSION(npart) :: nlrf_sph
double precision, intent(out), DIMENSION(npart) :: sqg

subroutine compute_pressure(npart, x, y, z, nlrf, eqos, pressure, verbose)

Wrapper function to compute the pressure from the given input

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)

Pressure at

logical, intent(in), optional :: verbose

subroutine correct_center_of_mass_of_system(npart, pos, nu, com_system)

Set the of the system to com_system

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npart
double precision, intent(inout) :: pos(3,npart)
double precision, intent(inout) :: nu(npart)
double precision, intent(in) :: com_system(3)

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

Wrapper function to compute the relativistic baryon mass density

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npart
double precision, intent(in) :: x(npart)
double precision, intent(in) :: y(npart)
double precision, intent(in) :: z(npart)
double precision, intent(in) :: nstar_sph(npart)
double precision, intent(out) :: nstar_id(npart)
double precision, intent(out) :: nlrf_sph(npart)
double precision, intent(out) :: sqg(npart)

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

Wrapper function to read the ID necessary to compute the relativistic baryonic mass

Arguments

Type IntentOptional Attributes Name
double precision, intent(in) :: x
double precision, intent(in) :: y
double precision, intent(in) :: z
double precision, intent(out) :: sqdetg
double precision, intent(out) :: baryon_density
double precision, intent(out) :: gamma_euler

subroutine integrate_mass_density(center, radius, central_density, dr, dth, dphi, mass, mass_profile, mass_profile_idx, radii, surf)

Wrapper function to integrate the relativistic baryonic mass density

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 place_particles_on_ellipsoidal_surfaces()

Place particles on ellipsoidal surfaces, and reflect particles with respect to the yz plane in the case of equal-mass binaries

FT 24.10.2022


Arguments

None

subroutine place_particles_on_lattices()

Place particles on lattices, one per matter object, and reflect particles with respect to the yz plane in the case of equal-mass binaries

FT 24.10.2022


Arguments

None

subroutine read_particles_from_formatted_file()

Read particles from formatted file, and reflect particles with respect to the yz plane in the case of equal-mass binaries

FT 21.10.2022


Arguments

None

subroutine read_particles_options()

Read the parameters in the file sphincs_id_particles.dat

FT 2022


Arguments

None

subroutine reflect_particles_yz_plane(pos_star1, pvol_star1, nu_star1, h_star1, npart_star1, pos_star2, pvol_star2, nu_star2, h_star2, npart_star2)

Reflect particles of star 1 with respect to the plane and place them on star 2

FT 07.02.2022


Arguments

Type IntentOptional Attributes Name
double precision, intent(in), DIMENSION(:,:) :: pos_star1

Array where to store the particle positions for star 1

double precision, intent(in), DIMENSION(:) :: pvol_star1

Array where to store the particle volumes for star 1

double precision, intent(in), DIMENSION(:) :: nu_star1

Array where to store the particle baryon number for star 1

double precision, intent(in), DIMENSION(:) :: h_star1

Array where to store the particle smoothing lengths for star 1

integer, intent(in) :: npart_star1

Variable where to store the particle number for star 1

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

Array where to store the particle positions for star 2

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

Array where to store the particle volumes for star 2

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

Array where to store the particle baryon number for star 2

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

Array where to store the particle smoothing lengths for star 2

integer, intent(inout) :: npart_star2

Variable where to store the particle number for star 2