perform_apm Module Procedure

module procedure perform_apm 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)

Uses

    • gradient
    • set_h
    • utility
    • units
    • matrix
    • sph_variables
    • APM
    • metric_on_particles
    • constants
    • RCB_tree_3D
    • analyze
    • sphincs_sph
    • kernel_table
  • proc~~perform_apm~~UsesGraph proc~perform_apm perform_apm APM APM proc~perform_apm->APM RCB_tree_3D RCB_tree_3D proc~perform_apm->RCB_tree_3D analyze analyze proc~perform_apm->analyze constants constants proc~perform_apm->constants gradient gradient proc~perform_apm->gradient kernel_table kernel_table proc~perform_apm->kernel_table matrix matrix proc~perform_apm->matrix metric_on_particles metric_on_particles proc~perform_apm->metric_on_particles module~utility utility proc~perform_apm->module~utility set_h set_h proc~perform_apm->set_h sph_variables sph_variables proc~perform_apm->sph_variables sphincs_sph sphincs_sph proc~perform_apm->sphincs_sph units units proc~perform_apm->units module~utility->constants module~utility->matrix

Compute the particle positions as follows:

  1. Take initial particle distribution as input
  2. Assume that the particles have the same mass
  3. Do the APM iteration so that the final SPH kernel estimate of the baryon mass density matches the baryon density in the given \(\mathrm{ID}\)
  4. Correct the particle masses ONCE, in order to match the density even better. Since we don't want a large mass ratio, we impose a maximum mass ratio when performing this correction.

After this procedure, the resulting particle distribution has positions and baryon numbers that kernel-estimate very well the mass density of the star, and has a low mass ratio.

This procedure assigns positions, smoothing lengths , and .

FT 04.06.2021


$OMP REDUCTION( MAX: radius_part ) &

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


Calls

proc~~perform_apm~~CallsGraph proc~perform_apm perform_apm 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 com com proc~perform_apm->com 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~check_particle_positions check_particle_positions proc~perform_apm->interface~check_particle_positions interface~correct_center_of_mass correct_center_of_mass proc~perform_apm->interface~correct_center_of_mass interface~impose_equatorial_plane_symmetry impose_equatorial_plane_symmetry proc~perform_apm->interface~impose_equatorial_plane_symmetry 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~is_finite_number is_finite_number proc~perform_apm->proc~is_finite_number proc~spherical_from_cartesian spherical_from_cartesian proc~perform_apm->proc~spherical_from_cartesian timer timer proc~perform_apm->timer proc~check_particle_positions check_particle_positions interface~check_particle_positions->proc~check_particle_positions proc~correct_center_of_mass correct_center_of_mass interface~correct_center_of_mass->proc~correct_center_of_mass proc~impose_equatorial_plane_symmetry impose_equatorial_plane_symmetry interface~impose_equatorial_plane_symmetry->proc~impose_equatorial_plane_symmetry 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->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 none~place_and_print_ghost_particles->proc~is_finite_number none~place_and_print_ghost_particles->proc~spherical_from_cartesian none~place_and_print_ghost_particles->timer 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 indexx indexx proc~check_particle_positions->indexx proc~correct_center_of_mass->com proc~correct_center_of_mass->proc~is_finite_number 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

Called by

proc~~perform_apm~~CalledByGraph proc~perform_apm perform_apm interface~perform_apm perform_apm interface~perform_apm->proc~perform_apm proc~construct_particles_std construct_particles_std proc~construct_particles_std->interface~perform_apm 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
integer, private :: a_numax
integer, private :: a_numax2
integer, private :: a_numin
integer, private :: a_numin2
double precision, private, DIMENSION(:,:), ALLOCATABLE :: all_pos
double precision, private, DIMENSION(:,:), ALLOCATABLE :: all_pos_prev
double precision, private, DIMENSION(:), ALLOCATABLE :: art_pr
double precision, private :: art_pr_max
double precision, private :: atmosphere_density
double precision, private :: average_displacement
integer, private :: b
integer, private :: cnt1
double precision, private, DIMENSION(:), ALLOCATABLE :: cnt_array
integer, private, DIMENSION(:), ALLOCATABLE :: cnt_move
integer, private :: cnt_push_ghost
double precision, private :: com_d
double precision, private :: com_x
double precision, private :: com_y
double precision, private :: com_z
double precision, private, DIMENSION(:,:), ALLOCATABLE :: correction_pos
double precision, private :: dN
double precision, private :: dN_av
double precision, private :: dN_max
double precision, private, DIMENSION(:), ALLOCATABLE :: dNstar
double precision, private, DIMENSION(:,:), ALLOCATABLE :: dS
double precision, private :: dS_norm_av
logical, private, parameter :: debug = .FALSE.
double precision, private :: dens_min
integer, private :: dim_seed
double precision, private :: dx
double precision, private :: dy
double precision, private :: dz
double precision, private :: ellipse_thickness
double precision, private :: err_N_max
double precision, private :: err_N_mean
double precision, private :: err_N_mean_min
double precision, private :: err_N_mean_min_old
double precision, private :: err_mean_old
double precision, private :: err_n_min
logical, private :: exist
character(len=:), private, ALLOCATABLE :: finalnamefile
type(timer), private :: find_h_bruteforce_timer
double precision, private :: frac
double precision, private :: ghost_displacement
double precision, private, DIMENSION(:,:), ALLOCATABLE :: ghost_pos
double precision, private :: h_av
double precision, private, DIMENSION(:), ALLOCATABLE :: h_guess
double precision, private, DIMENSION(:), ALLOCATABLE :: h_guess_tmp
double precision, private :: h_max
double precision, private, DIMENSION(:), ALLOCATABLE :: h_tmp
integer, private :: i
integer, private :: i_shell
integer, private :: ill
integer, private :: inde
integer, private :: index1
integer, private :: itot
integer, private :: itr
integer, private :: itr2
integer, private :: j
integer, private :: k
integer, private :: l
double precision, private :: l2norm_displacement
integer, private, parameter :: m_max_it = 50
double precision, private, parameter :: max_art_pr_ghost = 1.0D+10
integer, private, parameter :: max_npart = 20D+6
double precision, private :: max_nu
double precision, private :: max_nu2
integer, private :: max_push_ghost
double precision, private :: max_radius
double precision, private :: max_z_real
double precision, private :: mean_nu
double precision, private :: min_nu
double precision, private :: min_nu2
double precision, private :: min_radius
integer, private :: n_inc
integer, private, DIMENSION(:), ALLOCATABLE :: n_neighbors
integer, private :: n_problematic_h
double precision, private, DIMENSION(:,:), ALLOCATABLE :: nearest_neighbors
integer, private, DIMENSION(:), ALLOCATABLE :: neighbors_lists
double precision, private, DIMENSION(:), ALLOCATABLE :: nlrf_id
double precision, private :: nlrf_id_av
double precision, private, DIMENSION(:), ALLOCATABLE :: nlrf_sph
integer, private, parameter :: nn_des = 301
integer, private :: npart_all
integer, private :: npart_ghost
integer, private :: npart_real
integer, private :: npart_real_half
double precision, private, DIMENSION(:), ALLOCATABLE :: nstar_id
double precision, private :: nstar_id_av
double precision, private :: nstar_id_err
double precision, private, DIMENSION(:), ALLOCATABLE :: nstar_int
double precision, private, DIMENSION(:), ALLOCATABLE :: nstar_sph
double precision, private :: nstar_sph_av
double precision, private :: nstar_sph_err
double precision, private :: nu_all
double precision, private :: nu_ghost
double precision, private, DIMENSION(:), ALLOCATABLE :: nu_one
double precision, private :: nu_ratio
double precision, private, DIMENSION(:), ALLOCATABLE :: nu_tmp
double precision, private :: nu_tmp2
double precision, private :: nu_tot
integer, private :: nuratio_cnt
integer, private, parameter :: nuratio_max_steps = 300
integer, private, parameter :: nuratio_min_it = 100
double precision, private :: nuratio_tmp
double precision, private :: nuratio_tmp_prev
double precision, private, parameter :: nuratio_tol = 2.5D-3
integer, private :: nx
integer, private :: ny
integer, private :: nz
double precision, private, DIMENSION(:), ALLOCATABLE :: particle_density_final
double precision, private :: phi
double precision, private :: phi_ell
double precision, private, DIMENSION(:,:), ALLOCATABLE :: pos
double precision, private, DIMENSION(3) :: pos_corr_tmp
double precision, private, DIMENSION(3) :: pos_maxerr
double precision, private, DIMENSION(:,:), ALLOCATABLE :: pos_tmp
double precision, private, DIMENSION(:), ALLOCATABLE :: pressure_id
double precision, private :: pressure_id_av
double precision, private, DIMENSION(:), ALLOCATABLE :: pressure_sph
logical, private :: push_away_ghosts
double precision, private, DIMENSION(:), ALLOCATABLE :: pvol_tmp
double precision, private :: r
double precision, private :: r_ell
double precision, private :: rad_x
double precision, private :: rad_y
double precision, private :: rad_z
double precision, private :: radius_part
double precision, private :: radius_x
double precision, private :: radius_y
double precision, private :: radius_z
double precision, private :: rand_num
double precision, private :: rand_num2
integer, private :: rel_sign
double precision, private, DIMENSION(:), ALLOCATABLE :: rho_tmp
integer, private, parameter :: search_pos = 10
integer, private, DIMENSION(:), ALLOCATABLE :: seed
double precision, private, DIMENSION(:), ALLOCATABLE :: sqg
double precision, private :: stddev_dN
double precision, private :: stddev_nu
double precision, private :: theta
double precision, private :: theta_ell
double precision, private, parameter :: tiny_real = 1.0D-10
double precision, private :: tmp
double precision, private, parameter :: tol = 1.0D-3
double precision, private :: variance_dN
double precision, private :: variance_nu
double precision, private :: x_ell
double precision, private :: xmax
double precision, private :: xmin
double precision, private :: xtemp
double precision, private :: y_ell
double precision, private :: ymax
double precision, private :: ymin
double precision, private :: ytemp
double precision, private :: z_ell
double precision, private :: zmax
double precision, private :: zmin
double precision, private :: ztemp

Functions

function validate_position_final(x, y, z) result(answer)

Returns validate_position( x, y, z ) if the latter is present, 0 otherwise

FT 22.09.2021


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 logical

validate_position( x, y, z ) if the latter is present, 0 otherwise


Subroutines

subroutine allocate_apm_fields(npart_real, npart_ghost)

Allocate the fields used during the APM iteration

FT 20.04.2022


Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npart_real
integer, intent(in) :: npart_ghost

subroutine compute_hydro_momentum()

Compute the \(\mathrm{SPH}\) canonical momentum

FT 17.06.2022


Arguments

None

subroutine discard_atmosphere(npart_real)

Remove the atmosphere

FT 20.04.2022


Arguments

Type IntentOptional Attributes Name
integer, intent(inout) :: npart_real

subroutine dump_apm_pos()

Print the positions of the real and ghost particles to a formatted file

FT 20.04.2022


Arguments

None

subroutine get_nstar_id_atm(npart_real, x, y, z, nstar_sph, nstar_id, nlrf_sph, sqg, use_atmosphere)

Return various densities and the determinant of the metric, needed in the APM iteration

FT 5.12.2021


Arguments

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

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

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

Array of coordinates

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

Array of coordinates

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

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 sacetime metric DOUBLE PRECISION, INTENT(OUT):: nstar_eul_id(npart_real) Array to store the computed proper baryon number density seen by the Eulerian observer

logical, intent(in) :: use_atmosphere

.TRUE. if an atmosphere should be used during the APM, to allow the real aprticles more freedom to move around and adjust; .FALSE. otherwise

subroutine place_and_print_ghost_particles()

Place ghost particles around the matter object, and print their positions together with the positions of the real particles to a formatted file

FT 20.04.2022


Arguments

None

subroutine read_pressure_id(npart_real, x, y, z, pressure_id)

Return the pressure providd by the \(\mathrm{ID}\) (not the one computed from the \(\mathrm{SPH}\) density), needed in the APM iteration

FT 6.12.2022


Arguments

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

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

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

Array of coordinates

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

Array of coordinates

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

Array of coordinates

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

Array to store the pressure read from the \(\mathrm{ID}\)

subroutine reallocate_output_fields(npart_real)

Reallocate the fields to be returned by perform_apm

FT 20.04.2022


Arguments

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