DART-Ray
Data Types | Functions/Subroutines | Variables
rt_routines Module Reference

Contains the main radiation transfer routines. In particular, the routines to calculate the ray-cell intersections, while a ray propagates through the 3D grid, and the corresponding contributions to the radiation field energy density and scattered radiation. More...

Data Types

type  list_en_sca
 Contains the single values of en_sca which have to be sent to other nodes. More...
 

Functions/Subroutines

subroutine prepare_p_src ()
 Finds the host cells for the stellar point sources. It sets the values for src_cell() and cell_src(). More...
 
subroutine fix_ccoord_p_src (pos, rel, cellsize)
 Shifts slightly the source position within the host cell in case it is located on the border of the host cell. More...
 
subroutine calc_total_luminosity
 Calculates the total stellar luminosity within the RT model and store it in lumcell(). It also sets the arrays tot_rad_en() and tot_rad_en_or(). More...
 
subroutine calc_total_luminosity_sca
 Calculates tot_rad_en() and lumcell() at the beginning of each scattering iteration. More...
 
subroutine rt_loop
 This is the main RT loop over the sources of radiation. It is used in all calculation phases where the radiation field is derived. More...
 
subroutine rt_loop_2d
 Contains the main RT loop over the sources of radiation for the RT 2D mode. The loop is perfomed twice because the first instance processes only the cells with x,y,z > 0 while the second processes those with x =0 OR y=0 OR z=0. Repeating the loops is necessary because only the off-axis cells provide contributions toe the radiation field that can be easily symmetrised. More...
 
subroutine assign_scaspe_temp_arr (i)
 Assigns to scaspe_temp_arr() the values needed during the current OpenMP thread loop chunk. More...
 
logical function skip_cell_2d_rt_loop (i)
 Used to skip cell processing within rt_loop_2D. TRUE when the cell i does not have to be processed within 2D loop i2d. More...
 
subroutine wait_end_thread_loop
 Waits until all threads and all MPI processes have completed the source loop and sent out the remaining en_sca elements. More...
 
subroutine allocate_rt_loop_arrays
 Allocates arrays used in the rt_loop(), rt_loop_2d() and rt_loop_iobs() subroutines. More...
 
subroutine deallocate_rt_loop_arrays
 Deallocates arrays used in the rt_loop(), rt_loop_2d() and rt_loop_iobs() subroutines. More...
 
subroutine main_dir_loop (src_id, sid, src_lum)
 Contains the loop over the HEALPix directions used by rt_loop() and rt_loop_2D(). The minimum HEALPix resolution for the ray angular density corresponds to nside = 4. If necessary, the ray angular density is automatically increased or decreased along the ray path according to the input values for bm_par(), bm_par_max() and bm_par_sca(). More...
 
subroutine rt_loop_iobs
 Performs the ray-tracing loop to calculate the output specific intensity arrays i_obs(). More...
 
subroutine ray_tracing (nc, al, dl, clvl, src_lum, dplane_ray, prev_ray, isel_ray ,ray_type, ipix, cc_old)
 Performs the ray-tracing calculation for a single ray. The ray propagates through the model until it reach model edge or a blocking condition is met. The blocking condition depends on the RT calculation phase (see subroutine deposit()). More...
 
subroutine set_iq_list (ray_intensity, flag_beam)
 Sets the subscript list iq() of the ray_intensity() elements higher than zero and sets ray_intensity_dep(). More...
 
real(kind=real64) function, dimension(0:lnum-1) src_lum_value (id)
 Assigns the source luminosity within rt_loop() depending on the rt_type() and the cell/point source ID number. More...
 
subroutine calc_src_ccoord (id)
 Assigns the source coordinates to src_ccoord() depending on rt_type(). More...
 
subroutine calc_cproj (al, dl, cproj)
 Calculates the ray-direction projections cproj() on the coordinate axis. More...
 
subroutine find_cc_new2 (cc_old, cc_new, cproj, prev_ray, out, clvl)
 Finds the starting cell for a new ray (which is not necessarily the last cell crossed by the blocked ray from which the new ray is derived). More...
 
subroutine calc_outcube (rel_vec, cellsize, outcube)
 Calculates outcube, a xyz vector where an element is equal to 1 if the corresponding component of rel_vec is larger than the cell side divided by 2. This means that the selected point is outside the cell. More...
 
subroutine pre_calc_ads_arr
 Precalculates the ads_arr() arrays for low values of nside() (<=256) Ads_arr() elements for larger nside values are calculated when needed during the RT loop. More...
 
subroutine calc_ffn_arr (al, dl, ipix)
 Calculates the array ffn_arr() needed to derive the angular distribution of the scattered light after a ray-cell intersection. More...
 
subroutine calc_ads_arr (sin_dl, cos_dl, sin_al, cos_al)
 Calculates the array ads_arr(), which contains the values of the factors cos(theta) within the Henyey-Greenstein function for a single ray direction and for all directions of the scaspe() arrays. More...
 
subroutine calc_psel (incvec, psel, pabs, isel, clvl, cc, cproj)
 Calculates the length of the current ray-cell intersection. More...
 
subroutine deposit (cc, al, dl, ipix, ray_intensity, length, nc, psel, ray_status, flag_beam)
 Adds the contributions to the radiation field energy density array and to the scattered luminosity array of the cell intersected by a ray. It also updates the specific intensity of the ray which is modified by the ray-cell intersection if the cell is dusty. Note that the ray_intensity() array in this subroutine is different from the one in the caller subroutine. More...
 
real function calc_dist_cell (nc, cc)
 Calculates the distance between the centres of two cells. More...
 
subroutine check_wall_hit (cc, wall_hit)
 Checks whather the cell cc is located beyond the walls set in the input (see e.g. x_wall_on() and x_wall_coord() ). Note that the wall is not effective during the point source ray-tracing loop. More...
 
subroutine set_walls
 Sets the values of the wall coordinates in model units. More...
 
subroutine process_scatt_rad (cc, al, dl, ipix, en_sca)
 Adds the contribution of the scattered light luminosity to the crossed cell scaspe_arr() and scaspe_tot_arr() arrays. Note that this operation is normally performed through array vectorization since it is much faster than the corresponding do loop. However, when running the code in parallel, there is the possibility of race condition when two rays cross the same cell simultaneously. In order to avoid this, a do loop is performed with $OMP ATOMIC when lock_cell() = 1. Note also that if for the elements of scaspe_arr() not stored locally, the en_sca elements are stored in en_sca_list until transfering to the corresponding node. More...
 
subroutine set_wavelength_index (il, k)
 Sets the right wavelength subscript depending whether no_communications() is TRUE or FALSE. More...
 
subroutine define_next_level (flag_k)
 Determines the next value of nside() in the next iteration for main_dir_loop(). More...
 
subroutine calc_max_npix
 Finds maximum value in npix_arr within the range of wavelength used in the current RT algorithm. More...
 
subroutine create_scaspe
 Allocates and initialises the scattered luminosity array scaspe() as well as the sin/cos arrays ( theta_sca(), phi_sca(), sin_theta_sca(), cos_theta_sca(), sin_phi_sca(), cos_phi_sca()) for the fixed HEALPix directions and observer directions considered in the scaspe() arrays. Since these values depend only on kp_sca at each wavelength and on the observer lines-of-sight, the arrays are not stored for each wavelength but only for each kp_sca value present in kp_sca_arr(). More...
 
subroutine prepare_scaspe_splitting
 Prepares the splitting of the scaspe arrays within different node memory. It affects only calculations using multiple cluster nodes. It calculates the number of wavelengths lnum_node() for the scaspe array stored locally. More...
 
subroutine create_i_obs
 Allocates and initialises specific intensity arrays i_obs() and i_obs_in(). More...
 
subroutine create_scaspe_tot
 Allocates and initialises the scaspe_tot() array. Unlike the scaspe() array, the scaspe_tot() includes all of the scattered light at each position. If printed on an output file, it can be later used to calculate scattered light maps at arbitrary line-of-sight directions. More...
 
subroutine assign_src_lum (theta, phi, src_lum)
 Assigns the src_lum() values for rays carrying scattered light. This is done by simply locating the right element of the scaspe array. It is called in three cases: 1) for rays used in the radiation field calculation (rt_type = 'scatt'), 2) for rays used to derive i_obs() for an external observer in the rt_algorithm 'iobs', and 3) for rays used to derive the i_obs_in() value for internal observers. This version contains an interpolatation scheme over the sphere, but it is not used for the moment. More...
 
subroutine set_en_lim
 Sets the en_lim ( $ f_U $) parameter, which checks if a ray carries negligible luminosity and can be blocked without significantly reducing the calculation accuracy. More...
 
subroutine find_theta_phi_obs_in (ro, rc, theta, phi)
 Finds the theta and phi angles for the direction from a source (emitting cell or point source) to an internal observer. More...
 
subroutine create_psel_av_arr
 Creates the arrays that stores the average path of the rays originating from each source. The array is expanded at the beginning of each scattering iteration/dust heating iteration in expand_psel_av_arr(). More...
 
subroutine expand_psel_av_arr
 Expands the psel_av_arr() array at the beginning of each scattering iteration and of each dust heating iteration. More...
 
subroutine check_grid_symmetry
 Checks that the 3D grid is indeed axisymmetric in the rt_algorithm 2D. More...
 
subroutine find_linked_cells (nc, link_list)
 Finds the linked cells in the 2D mode (that is, the cells located at the symmetric points). More...
 
subroutine check_2d_src
 Checks that no more than one source is allowed at the grid origin in the 2D mode. This is needed because it is the only case where symmetries can be exploited in the same way as for the emitting cells. More...
 
subroutine fix_symmetry
 Takes care of symmetrising the radiation fields and scaspe() arrays at the end of the RT loops in the RT 2D mode. R and z symmetry is exploited here. More...
 
subroutine fix_symmetry_part1
 Handles the array symmetrization after the RT precalculation in 2D mode. More...
 
subroutine fix_symmetry_part2
 Takes care of symmetrizing radiation fields and scaspe() arrays at the end of direct light processing for the RT 2D mode. More...
 
subroutine fix_symmetry_part3
 Symmetrises the radiation fields and scaspe() arrays at the end of each scattering iteration in the RT 2D mode. More...
 
subroutine remove_negative (array)
 Assigns zero to all negative elements of an array. Used in fix_symmetry_part3() to remove negative due to numerical accuracy. Remember that this could give problems in case you want to work with negative intensities. More...
 
subroutine remove_negative_2d (array)
 Same as remove_negative() but for 2D arrays. More...
 
subroutine create_scaspe_prev
 Allocates and initialises a new scaspe array ( scaspe_prev_arr()) which is used while symmetrising the scaspe_arr() array in the RT 2D mode or in the sequential scattering algorithm (see sequential_scattering()). More...
 
subroutine calc_scaspe_indices
 Derives the re-ordering scaspe() indices for symmetrising the scaspe arrays in the RT 2D mode. More...
 
subroutine create_en_sca_list
 Allocates and initializes en_sca_list() array and the counter count_en_sca(). More...
 
subroutine add_to_en_sca_list (cc, ipix, il, en_sca_out)
 Adds additional element to the en_sca_list() of the OPENMP thread. If the list contains the maximum number of elements allowed, the transfering subroutine handle_mpi_transfer is called. More...
 
subroutine process_en_sca_list
 Processes the en_sca_list arrays of each OPENMP thread. It joins them in a single array. Then the master thread of each MPI process handles the communication and processing. More...
 
subroutine sort_en_sca_list (el_out, i0_el_out)
 Sorts the elements of en_sca_list for each thread before passing it to en_sca_list_all(). This avoids to make the sorting later within the master block. More...
 
subroutine process_en_sca_received (im)
 Processes en_sca_list_received within each MPI process. This is done by all threads bu the master thread. More...
 
subroutine prepare_en_sca_transfer
 Prepares en_sca_list_all() and determines how many elements to send to each MPI process and starting indeces. More...
 
subroutine en_sca_transfer ()
 Transfers en_sca_list elements between the MPI processes and process them. The transfering is done by the master threads while the other threads process the elements as soon as they arrive. This method improves efficiency because communication and computation are partly overlapped. More...
 
subroutine create_mpi_type
 Creates the MPI data type used for transfering the en_sca_list elements between MPI processes. More...
 
subroutine assign_scaspe_temp (cc)
 Assigns values to scaspe_temp_arr() array. It calls the subroutines to send and receive scaspe_arr() arrays values from the other MPI processes. More...
 
subroutine prepare_scaspe_temp_transfer
 Prepares the scaspe_temp_send() array whose parts will be transfered to other nodes using mpi_scatterv. To reduce the number of scaspe transfers, a group of scaspe_temp_arr arrays (corresponding to many cells not just those requested) are passed to the MPI process/OPENMP Thread that sent the request. This group corresponds to the scaspe_temp_arr elements of the cells that will be processed sequentially by that MPI process/OPENMP thread within the current OPENMP chunk (see OMP DO SCHEDULE). The order in which elements are assigned to the 1-dim scaspe_temp_send() array is angular directions (npix_arr(kl)), wavelength(kl), cell id (ic) More...
 
subroutine scaspe_temp_transfer
 Sends chunks of scaspe_temp_send() to the other processes using mpi_scatterv. It stores the received values in scaspe_temp_recv(). More...
 
subroutine handle_mpi_transfers (transfer_type, cc)
 Handles MPI transfers depending on the transfer_type input from each OPENMP thread and MPI PROCESS. The reason for this subroutine is that there are two different situations where a thread can ask for communication to the other MPI process: 1) the assignment of scaspe_temp_arr() arrays; 2) the transfer and processing of the en_sca_list() elements. This routine is called when all the OPENMP threads within an MPI process have reached one of the two calling points (depending on the operation needed). When all the MPI processes have arrived to mpi_allreduce, they communicate transfer_type_tot_local(), which is equal to the number of threads only if no thread needs communication of scaspe_temp elements (that is, all threads requested communication of en_sca_list elements). Then the routines assign_scaspe_temp() and process_en_sca_list() are called to perform the communications. More...
 
subroutine count_processed_cells
 It counts the number of cells that have been processed so far. Only those with non zero luminosity are counted. This count can be used when estimating the average processing time for a single cell. More...
 
subroutine check_im (im)
 It checks whether the thread is processing a cell with an im value which is the start of an OPENMP loop chunk. If so, it sets iscaspe_big = num_scaspe_pass so that the query of scaspe_temp_arr_big is made in the rt_loop. Otherwise, it just increases the value of iscaspe_big. This subroutine is not important in the no communication mode (see no_communications()). More...
 
subroutine scale_dens_arr
 Assigns the values of dens_arr() by scaling the dens() array read from the main grid file. More...
 
subroutine set_i_opacity_arrays (i0, i1)
 Sets the subscript of the first element of the value segment within kext_arr(), gsca_arr() and ksca_norm_arr() used in the RT calculation. More...
 
subroutine put_lock_to_cell (i)
 Increases lock_cell() by a unity. More...
 
subroutine remove_lock_from_cell (i)
 Decreases lock_cell() by a unity. More...
 
subroutine store_reshape_arrays
 Changes the size of the wavelength dimension of u_fest_arr() to match the dust emission wavelength range. Stores the radiation field energy density in the stellar emission wavelength range into u_final_uv_opt(). New u_final_arr and dens_arr arrays are created to cover the wavelength range of the dust emission. Note that in this subroutine, lnum is still set to lnum_stars. More...
 
subroutine set_units
 Sets units_i_obs(), units_ufield() and cs(). More...
 
subroutine update_i_obs_arr
 Adds the i_obs_arr(), calculated durint the current dust heating iteration, to the arrays i_obs_arr_dir() and i_obs_arr_tot(). More...
 
subroutine set_npix_arr
 Sets the number of pixels for the scaspe() arrays at each wavelength (npix_arr() and npix_hp_arr()). Having a variable size for the scattering source function allows to save memory especially in the infrared. More...
 
subroutine assign_i_obs_to_project
 Assigns i_obs_arr() values for the projection algorithm. More...
 

Variables

integer(kind=int32) kp_sca_max
 
integer(kind=int32) nside_sca
 
integer(kind=int32) bm_par
 
integer(kind=int32) bm_par_sca
 
integer(kind=int32) bm_par_max
 
integer(kind=int32) pabs_max
 
real(kind=real64) rad_lim
 
real(kind=real64) en_lim
 
real(kind=real64) conv_en_lim
 
real(kind=real64) tau_cell_max
 
real(kind=real64) accuracy
 
real(kind=real64), dimension(:), allocatable tot_rad_en
 
real(kind=real64), dimension(:), allocatable tot_rad_en_or
 
integer rt_type
 
character(len=20) rt_algorithm
 
integer rt_algorithm_id
 
logical done
 
logical cnflag
 
logical cnflag_dust
 
integer nside
 
integer clvl_old
 
integer, parameter nside_min = 4
 
integer, parameter npix_main =12
 
real(kind=real64), dimension(3) src_ccoord
 
integer, dimension(:,:), allocatable ccindd_nc
 
integer, dimension(:,:), allocatable ccindd
 
integer nproc
 
integer ipsel_av
 
integer ipsel_av_tot
 
real(kind=real64), parameter glepsilon =1.0e-7
 
type(var_arr_2d), dimension(:), allocatable ads_arr4
 
type(var_arr_2d), dimension(:), allocatable ads_arr8
 
type(var_arr_2d), dimension(:), allocatable ads_arr16
 
type(var_arr_2d), dimension(:), allocatable ads_arr32
 
type(var_arr_2d), dimension(:), allocatable ads_arr64
 
type(var_arr_2d), dimension(:), allocatable ads_arr128
 
type(var_arr_2d), dimension(:), allocatable ads_arr256
 
type(var_arr_1d), dimension(:), allocatable theta_sca
 
type(var_arr_1d), dimension(:), allocatable phi_sca
 
type(var_arr_1d), dimension(:), allocatable sin_theta_sca
 
type(var_arr_1d), dimension(:), allocatable sin_phi_sca
 
type(var_arr_1d), dimension(:), allocatable cos_theta_sca
 
type(var_arr_1d), dimension(:), allocatable cos_phi_sca
 
type(var_arr_1d), dimension(:), allocatable ffn_arr_mw
 
type(var_arr_1d), dimension(:), allocatable ffn_arr
 
type(var_arr_1d), dimension(:), allocatable ads_arr
 
real(kind=real64) cs
 
real(kind=real64) psel_av
 
real(kind=real64), dimension(3) pabs_arr
 
real(kind=real64), dimension(:), allocatable i_obs_temp
 
real(kind=real64), dimension(:), allocatable lum_lost_temp
 
real(kind=real64), dimension(:), allocatable lum_lost
 
real(kind=real64), dimension(:), allocatable lum_lost_prev
 
integer iterations
 
integer iterations_dustem
 
real(kind=real64) vec_mod
 
type(var_arr_1d), dimension(:), allocatable ix
 
type(var_arr_1d), dimension(:), allocatable iy
 
type(var_arr_1d), dimension(:), allocatable iz
 
type(var_arr_1d), dimension(:), allocatable ixy
 
type(var_arr_1d), dimension(:), allocatable ixz
 
type(var_arr_1d), dimension(:), allocatable iyz
 
type(var_arr_1d), dimension(:), allocatable ixyz
 
integer chunk_size
 
integer, dimension(:), allocatable npix_arr
 
integer npix_hp
 
integer, dimension(:), allocatable tot_npix_arr_local
 
integer, dimension(:), allocatable npix_hp_arr
 
integer, dimension(:), allocatable npix_unique
 
integer, dimension(:), allocatable npix_hp_unique
 
integer, dimension(:), allocatable kp_unique
 
integer dim_npix_unique
 
integer max_npix
 
integer, dimension(:), allocatable kp_sca_arr
 
integer, dimension(:), allocatable ik_sca_arr
 
integer i2d
 
integer, parameter rtt_start = 0
 
integer, parameter rtt_precalc_cell = 1
 
integer, parameter rtt_precalc_src = 2
 
integer, parameter rtt_output_part1 = 3
 
integer, parameter rtt_prep_part2 = 4
 
integer, parameter rtt_dir_cell = 5
 
integer, parameter rtt_i_obs_dir_cell = 6
 
integer, parameter rtt_dir_src = 7
 
integer, parameter rtt_i_obs_dir_src = 8
 
integer, parameter rtt_output_part2 = 9
 
integer, parameter rtt_scatt = 10
 
integer, parameter rtt_i_obs = 11
 
integer, parameter rtt_read_i_obs_part2 = 12
 
integer, parameter rtt_read_i_obs = 13
 
integer, parameter rtt_read_ufield = 14
 
integer, parameter rtt_read_scaspe_tot = 15
 
integer, parameter rtt_grid_init_stars = 16
 
integer, parameter rtt_grid_init_dust = 17
 
integer, parameter rtt_grid_init_projection = 18
 
integer, parameter rta_main = 0
 
integer, parameter rta_i_obs = 1
 
integer, parameter rta_2d = 2
 
integer, parameter rta_sed = 3
 
integer, parameter rta_dust = 4
 
integer, parameter rta_dust2d = 5
 
integer, parameter rta_i_obs_dust = 6
 
integer, parameter rta_sed_dust = 7
 
integer, parameter rta_projection = 8
 
integer, parameter ras_go_high = 0
 
integer, parameter ras_gone = 1
 
integer, parameter ras_go_low = 2
 
integer, parameter ras_first_launch = 3
 
integer, parameter ras_re_launched = 4
 
logical print_scaspe_tot
 
logical print_psel_av
 
logical sequential_scattering
 
real(kind=real64), parameter tol_p = 1E-6
 param tol_p_src Tolerance parameter to find host cell of a point source More...
 
logical, dimension(:), allocatable iq_a
 
integer, dimension(:), allocatable iq
 
integer lnum_a
 
integer lnum_a_old
 
integer lnum_node
 
integer lnum_maps
 
integer lnum_node_maps
 
integer, dimension(:), allocatable lnum_node_arr
 
logical, dimension(:), allocatable iq_sca_node
 
integer, dimension(:), allocatable iq_sca_id
 
integer, dimension(:), allocatable iq_maps_id
 
integer, dimension(:), allocatable im_lambda_arr
 
real(kind=real64), dimension(:), allocatable ray_intensity_dep
 
logical no_communications
 
type(list_en_sca), dimension(:), allocatable en_sca_list_received
 
integer, dimension(:), allocatable en_sca_id_thread
 
type(list_en_sca), dimension(:), allocatable en_sca_list
 
type(list_en_sca), dimension(:), allocatable en_sca_list_all
 
type(list_en_sca), dimension(:), allocatable temp_en_sca_list
 
integer en_sca_arrtype
 
integer, dimension(:), allocatable ind_en_sca_list
 
integer count_en_sca
 
integer, dimension(:), allocatable count_en_sca_arr
 
integer count_en_sca_tot
 
integer, dimension(:), allocatable i0_count_en_sca
 
integer, dimension(:,:), allocatable el_out_arr
 
integer, dimension(:,:), allocatable i0_el_out_arr
 
integer, dimension(:), allocatable el_out_mpi
 
integer, dimension(:), allocatable el_in_mpi
 
integer, dimension(:), allocatable i0_el_out_mpi
 
integer, dimension(:), allocatable i0_el_in_mpi
 
logical, dimension(:), allocatable en_sca_confirm
 
integer size_en_sca_list
 
integer n_el_in_tot
 
integer threadid
 
logical, dimension(:), allocatable mpi_proc_completed
 
logical, dimension(:), allocatable wait_thread_arr
 
logical, dimension(:), allocatable handling_mpi_thread_arr
 
integer count_transfer
 
integer, parameter count_transfer_start =100
 
integer num_processed_cells
 
integer, dimension(:,:), allocatable cc_list_all
 
integer, dimension(:), allocatable cc_list_local
 
type(var_arr_2d), dimension(:), allocatable scaspe_temp_recv
 
real(kind=real64), dimension(:), allocatable scaspe_temp_send
 
type(var_arr_1d), dimension(:), allocatable scaspe_temp_arr
 
type(var_arr_2d), dimension(:), allocatable scaspe_temp_arr_big
 
integer iscaspe_big
 
integer *1, dimension(:,:), allocatable transfer_type_all
 
integer *1, dimension(:), allocatable transfer_type_local
 
integer *1 transfer_type_tot_local
 
integer transfer_type_tot_all
 
integer num_scaspe_pass
 
integer *1, parameter mpi_transfer_scaspe = 0
 
integer *1, parameter mpi_transfer_en_sca = 1
 
real(kind=real64) dist_obs
 
integer, dimension(:), allocatable ind_i_obs
 
integer, dimension(:), allocatable ind_out_maps
 
integer, parameter max_n_input = 1000
 
real(kind=real64), dimension(:,:), allocatable sed_arr
 
real(kind=real64), dimension(:,:), allocatable sed_arr_dir
 
logical test_run
 
logical x_wall_on
 
logical y_wall_on
 
logical z_wall_on
 
real(kind=real64), dimension(2) x_wall_coord
 
real(kind=real64), dimension(2) y_wall_coord
 
real(kind=real64), dimension(2) z_wall_coord
 
logical only_direct_rt
 
character(len=lcar) param_to_project
 

Detailed Description

Contains the main radiation transfer routines. In particular, the routines to calculate the ray-cell intersections, while a ray propagates through the 3D grid, and the corresponding contributions to the radiation field energy density and scattered radiation.

Function/Subroutine Documentation

subroutine rt_routines::add_to_en_sca_list ( integer  cc,
integer  ipix,
integer  il,
real(kind=real64)  en_sca_out 
)

Adds additional element to the en_sca_list() of the OPENMP thread. If the list contains the maximum number of elements allowed, the transfering subroutine handle_mpi_transfer is called.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::allocate_rt_loop_arrays ( )

Allocates arrays used in the rt_loop(), rt_loop_2d() and rt_loop_iobs() subroutines.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::assign_i_obs_to_project ( )

Assigns i_obs_arr() values for the projection algorithm.

Here is the caller graph for this function:

subroutine rt_routines::assign_scaspe_temp ( integer  cc)

Assigns values to scaspe_temp_arr() array. It calls the subroutines to send and receive scaspe_arr() arrays values from the other MPI processes.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::assign_scaspe_temp_arr ( integer  i)

Assigns to scaspe_temp_arr() the values needed during the current OpenMP thread loop chunk.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::assign_src_lum ( real(kind=real64), intent(in)  theta,
real(kind=real64), intent(in)  phi,
real(kind=real64), dimension(0:lnum-1)  src_lum 
)

Assigns the src_lum() values for rays carrying scattered light. This is done by simply locating the right element of the scaspe array. It is called in three cases: 1) for rays used in the radiation field calculation (rt_type = 'scatt'), 2) for rays used to derive i_obs() for an external observer in the rt_algorithm 'iobs', and 3) for rays used to derive the i_obs_in() value for internal observers. This version contains an interpolatation scheme over the sphere, but it is not used for the moment.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::calc_ads_arr ( real(kind=real64)  sin_dl,
real(kind=real64)  cos_dl,
real(kind=real64)  sin_al,
real(kind=real64)  cos_al 
)

Calculates the array ads_arr(), which contains the values of the factors cos(theta) within the Henyey-Greenstein function for a single ray direction and for all directions of the scaspe() arrays.

Parameters
[in]sin_dlSine of the theta angle of the ray direction
[in]cos_dlCosine of the theta angle of the ray direction
[in]sin_alSine of the phi angle of the ray direction
[in]cos_alCosine of the phi angle of the ray direction

Here is the caller graph for this function:

subroutine rt_routines::calc_cproj ( real(kind=real64)  al,
real(kind=real64)  dl,
real(kind=real64), dimension(3)  cproj 
)

Calculates the ray-direction projections cproj() on the coordinate axis.

Parameters
[in]alangle in rad between ray-projection on XY plane and X axis
[in]dlangle in rad between ray and Z axis
[out]cprojunit vector with ray direction components

Here is the caller graph for this function:

real function rt_routines::calc_dist_cell ( integer  nc,
integer  cc 
)

Calculates the distance between the centres of two cells.

subroutine rt_routines::calc_ffn_arr ( real(kind=real64)  al,
real(kind=real64)  dl,
integer  ipix 
)

Calculates the array ffn_arr() needed to derive the angular distribution of the scattered light after a ray-cell intersection.

Parameters
[in]alPhi angle of the ray direction
[in]dlTheta angle of the ray direction
[in]ipixHEALPix number of the ray direction

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::calc_max_npix ( )

Finds maximum value in npix_arr within the range of wavelength used in the current RT algorithm.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::calc_outcube ( real(kind=real64), dimension(3)  rel_vec,
real(kind=real64)  cellsize,
integer, dimension(3)  outcube 
)

Calculates outcube, a xyz vector where an element is equal to 1 if the corresponding component of rel_vec is larger than the cell side divided by 2. This means that the selected point is outside the cell.

Parameters
[in]rel_vecVector from cell centre to considered point
[in]cellsizeSize of the cell
[out]outcubeVector whose element are equal to 1 if the projection of rel_vec on the corresponding axis is larger than cellsize

Here is the caller graph for this function:

subroutine rt_routines::calc_psel ( integer, dimension(3)  incvec,
real(kind=real64)  psel,
real(kind=real64), dimension(3)  pabs,
integer  isel,
integer  clvl,
integer  cc,
real(kind=real64), dimension(3)  cproj 
)

Calculates the length of the current ray-cell intersection.

Parameters
[in]incvecxyz vector equal to +1 or -1 depending if the ray is pointing towards the positive or negative direction of the corresponding axis. It is equal to 0 if the ray is along a perpendicular axis.
[out]pselLength of the ray-cell intersection
[out]pabsArray of possible lengths for the ray-cell intersection. These are the lengths in the cases the ray hits planes perpendicular to x,y or z directions when exiting from the intersected cell.
[out]iselDirection perpendicular to the plane hit by the ray when exiting the intersected cell (0 = x, 1= y, 2=z).
[in]clvlSubdivision level of the intersected cell
[in]ccID number of the intersected cell
[in]cprojRay direction projections on coordinate axis

Here is the caller graph for this function:

subroutine rt_routines::calc_scaspe_indices ( )

Derives the re-ordering scaspe() indices for symmetrising the scaspe arrays in the RT 2D mode.

Here is the caller graph for this function:

subroutine rt_routines::calc_src_ccoord ( integer  id)

Assigns the source coordinates to src_ccoord() depending on rt_type().

Here is the caller graph for this function:

subroutine rt_routines::calc_total_luminosity ( )

Calculates the total stellar luminosity within the RT model and store it in lumcell(). It also sets the arrays tot_rad_en() and tot_rad_en_or().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::calc_total_luminosity_sca ( )

Calculates tot_rad_en() and lumcell() at the beginning of each scattering iteration.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::check_2d_src ( )

Checks that no more than one source is allowed at the grid origin in the 2D mode. This is needed because it is the only case where symmetries can be exploited in the same way as for the emitting cells.

Here is the caller graph for this function:

subroutine rt_routines::check_grid_symmetry ( )

Checks that the 3D grid is indeed axisymmetric in the rt_algorithm 2D.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::check_im ( integer  im)

It checks whether the thread is processing a cell with an im value which is the start of an OPENMP loop chunk. If so, it sets iscaspe_big = num_scaspe_pass so that the query of scaspe_temp_arr_big is made in the rt_loop. Otherwise, it just increases the value of iscaspe_big. This subroutine is not important in the no communication mode (see no_communications()).

Here is the caller graph for this function:

subroutine rt_routines::check_wall_hit ( integer  cc,
logical  wall_hit 
)

Checks whather the cell cc is located beyond the walls set in the input (see e.g. x_wall_on() and x_wall_coord() ). Note that the wall is not effective during the point source ray-tracing loop.

Here is the caller graph for this function:

subroutine rt_routines::count_processed_cells ( )

It counts the number of cells that have been processed so far. Only those with non zero luminosity are counted. This count can be used when estimating the average processing time for a single cell.

Here is the caller graph for this function:

subroutine rt_routines::create_en_sca_list ( )

Allocates and initializes en_sca_list() array and the counter count_en_sca().

Here is the caller graph for this function:

subroutine rt_routines::create_i_obs ( )

Allocates and initialises specific intensity arrays i_obs() and i_obs_in().

Here is the caller graph for this function:

subroutine rt_routines::create_mpi_type ( )

Creates the MPI data type used for transfering the en_sca_list elements between MPI processes.

Here is the caller graph for this function:

subroutine rt_routines::create_psel_av_arr ( )

Creates the arrays that stores the average path of the rays originating from each source. The array is expanded at the beginning of each scattering iteration/dust heating iteration in expand_psel_av_arr().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::create_scaspe ( )

Allocates and initialises the scattered luminosity array scaspe() as well as the sin/cos arrays ( theta_sca(), phi_sca(), sin_theta_sca(), cos_theta_sca(), sin_phi_sca(), cos_phi_sca()) for the fixed HEALPix directions and observer directions considered in the scaspe() arrays. Since these values depend only on kp_sca at each wavelength and on the observer lines-of-sight, the arrays are not stored for each wavelength but only for each kp_sca value present in kp_sca_arr().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::create_scaspe_prev ( )

Allocates and initialises a new scaspe array ( scaspe_prev_arr()) which is used while symmetrising the scaspe_arr() array in the RT 2D mode or in the sequential scattering algorithm (see sequential_scattering()).

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::create_scaspe_tot ( )

Allocates and initialises the scaspe_tot() array. Unlike the scaspe() array, the scaspe_tot() includes all of the scattered light at each position. If printed on an output file, it can be later used to calculate scattered light maps at arbitrary line-of-sight directions.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::deallocate_rt_loop_arrays ( )

Deallocates arrays used in the rt_loop(), rt_loop_2d() and rt_loop_iobs() subroutines.

Here is the caller graph for this function:

subroutine rt_routines::define_next_level ( logical  flag_k)

Determines the next value of nside() in the next iteration for main_dir_loop().

Parameters
[out]flag_kLogical variable equal to TRUE if there are no rays left in the ray lists.

Here is the caller graph for this function:

subroutine rt_routines::deposit ( integer, intent(in)  cc,
real(kind=real64)  al,
real(kind=real64)  dl,
integer, intent(in)  ipix,
real(kind=real64), dimension(0:lnum_a-1)  ray_intensity,
real(kind=real64), intent(in)  length,
integer, intent(in)  nc,
real(kind=real64)  psel,
integer  ray_status,
logical  flag_beam 
)

Adds the contributions to the radiation field energy density array and to the scattered luminosity array of the cell intersected by a ray. It also updates the specific intensity of the ray which is modified by the ray-cell intersection if the cell is dusty. Note that the ray_intensity() array in this subroutine is different from the one in the caller subroutine.

Parameters
[in]ccID number of the intersected cell
[in,out]ray_intensitySpecific intensity of the ray
[in]lengthLength of the ray-cell intersection
[in]ncID number of the cell originating the ray
[in]pselPath crossed by the ray
[in,out]ray_statusLabel indicating whether the ray has to be blocked and if the angular density of the rays has to be increased or decreased. Possible values ('first_launch' ,'re-launched', 'gone', 'go_high' , 'go_low')
[in]flag_beamFlag with TRUE value when the ray has just started from the original cell or has been just relaunched from an intermediate point. It avoids infinite loops when new angular resolution is too low after being too high.!>
[in]ray_typeRay_type ID

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::en_sca_transfer ( )

Transfers en_sca_list elements between the MPI processes and process them. The transfering is done by the master threads while the other threads process the elements as soon as they arrive. This method improves efficiency because communication and computation are partly overlapped.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::expand_psel_av_arr ( )

Expands the psel_av_arr() array at the beginning of each scattering iteration and of each dust heating iteration.

Here is the caller graph for this function:

subroutine rt_routines::find_cc_new2 ( integer  cc_old,
integer  cc_new,
real(kind=real64), dimension(3)  cproj,
real(kind=real64)  prev_ray,
logical  out,
integer  clvl 
)

Finds the starting cell for a new ray (which is not necessarily the last cell crossed by the blocked ray from which the new ray is derived).

Parameters
[in]cc_oldID number of the last cell intersected by the blocked ray
[out]cc_newID number of the cell intersected by the new ray
[in]cprojRay direction projections on coordinate axis
[in]prev_rayDistance crossed by the blocked ray
[out]outFlag equal to TRUE if ray is now outside the model
[out]clvlSubdivision level of the cell intersected by the new ray

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::find_linked_cells ( integer  nc,
integer, dimension(0:6)  link_list 
)

Finds the linked cells in the 2D mode (that is, the cells located at the symmetric points).

Parameters
[in]ncID number of a cell
[out]link_listArray of ID number of the linked cells

Here is the caller graph for this function:

subroutine rt_routines::find_theta_phi_obs_in ( real(kind= real64), dimension(3)  ro,
real(kind= real64), dimension(3)  rc,
real(kind= real64)  theta,
real(kind= real64)  phi 
)

Finds the theta and phi angles for the direction from a source (emitting cell or point source) to an internal observer.

Parameters
[in]iID number of the internal observer
[out]thetaSource - observer direction theta angle
[out]phiSource - observer direction phi angle

Here is the caller graph for this function:

subroutine rt_routines::fix_ccoord_p_src ( real(kind=real64), dimension(3)  pos,
real(kind=real64), dimension(3)  rel,
real(kind=real64)  cellsize 
)

Shifts slightly the source position within the host cell in case it is located on the border of the host cell.

Here is the caller graph for this function:

subroutine rt_routines::fix_symmetry ( )

Takes care of symmetrising the radiation fields and scaspe() arrays at the end of the RT loops in the RT 2D mode. R and z symmetry is exploited here.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::fix_symmetry_part1 ( )

Handles the array symmetrization after the RT precalculation in 2D mode.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::fix_symmetry_part2 ( )

Takes care of symmetrizing radiation fields and scaspe() arrays at the end of direct light processing for the RT 2D mode.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::fix_symmetry_part3 ( )

Symmetrises the radiation fields and scaspe() arrays at the end of each scattering iteration in the RT 2D mode.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::handle_mpi_transfers ( integer*1  transfer_type,
integer  cc 
)

Handles MPI transfers depending on the transfer_type input from each OPENMP thread and MPI PROCESS. The reason for this subroutine is that there are two different situations where a thread can ask for communication to the other MPI process: 1) the assignment of scaspe_temp_arr() arrays; 2) the transfer and processing of the en_sca_list() elements. This routine is called when all the OPENMP threads within an MPI process have reached one of the two calling points (depending on the operation needed). When all the MPI processes have arrived to mpi_allreduce, they communicate transfer_type_tot_local(), which is equal to the number of threads only if no thread needs communication of scaspe_temp elements (that is, all threads requested communication of en_sca_list elements). Then the routines assign_scaspe_temp() and process_en_sca_list() are called to perform the communications.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::main_dir_loop ( integer  src_id,
integer  sid,
real(kind=real64), dimension(0:lnum-1)  src_lum 
)

Contains the loop over the HEALPix directions used by rt_loop() and rt_loop_2D(). The minimum HEALPix resolution for the ray angular density corresponds to nside = 4. If necessary, the ray angular density is automatically increased or decreased along the ray path according to the input values for bm_par(), bm_par_max() and bm_par_sca().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::pre_calc_ads_arr ( )

Precalculates the ads_arr() arrays for low values of nside() (<=256) Ads_arr() elements for larger nside values are calculated when needed during the RT loop.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::prepare_en_sca_transfer ( )

Prepares en_sca_list_all() and determines how many elements to send to each MPI process and starting indeces.

Here is the caller graph for this function:

subroutine rt_routines::prepare_p_src ( )

Finds the host cells for the stellar point sources. It sets the values for src_cell() and cell_src().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::prepare_scaspe_splitting ( )

Prepares the splitting of the scaspe arrays within different node memory. It affects only calculations using multiple cluster nodes. It calculates the number of wavelengths lnum_node() for the scaspe array stored locally.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::prepare_scaspe_temp_transfer ( )

Prepares the scaspe_temp_send() array whose parts will be transfered to other nodes using mpi_scatterv. To reduce the number of scaspe transfers, a group of scaspe_temp_arr arrays (corresponding to many cells not just those requested) are passed to the MPI process/OPENMP Thread that sent the request. This group corresponds to the scaspe_temp_arr elements of the cells that will be processed sequentially by that MPI process/OPENMP thread within the current OPENMP chunk (see OMP DO SCHEDULE). The order in which elements are assigned to the 1-dim scaspe_temp_send() array is angular directions (npix_arr(kl)), wavelength(kl), cell id (ic)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::process_en_sca_list ( )

Processes the en_sca_list arrays of each OPENMP thread. It joins them in a single array. Then the master thread of each MPI process handles the communication and processing.

Parameters
el_outNumber of en_sca elements to be sent to each MPI process.
i0_el_outList of positions of the first elements in the en_sca_list_all() array to be sent to each MPI process.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::process_en_sca_received ( integer  im)

Processes en_sca_list_received within each MPI process. This is done by all threads bu the master thread.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::process_scatt_rad ( integer  cc,
real (kind=real64)  al,
real (kind=real64)  dl,
integer  ipix,
real (kind=real64), dimension(0:lnum_a-1)  en_sca 
)

Adds the contribution of the scattered light luminosity to the crossed cell scaspe_arr() and scaspe_tot_arr() arrays. Note that this operation is normally performed through array vectorization since it is much faster than the corresponding do loop. However, when running the code in parallel, there is the possibility of race condition when two rays cross the same cell simultaneously. In order to avoid this, a do loop is performed with $OMP ATOMIC when lock_cell() = 1. Note also that if for the elements of scaspe_arr() not stored locally, the en_sca elements are stored in en_sca_list until transfering to the corresponding node.

Parameters
[in]ccID number of the cell intersected by a ray
[in]en_ScaTotal luminosity scattered by the intersected cell
[in]iwWavelength index

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::put_lock_to_cell ( integer  i)

Increases lock_cell() by a unity.

Here is the caller graph for this function:

subroutine rt_routines::ray_tracing ( integer  nc,
real(kind=real64), intent(in)  al,
real(kind=real64), intent(in)  dl,
integer  clvl,
real(kind=real64), dimension(0:lnum-1), intent(in)  src_lum,
real(kind=real64)  dplane_ray,
real(kind=real64)  prev_ray,
integer  isel_ray,
integer  ray_type,
integer  ipix,
integer  cc_old 
)

Performs the ray-tracing calculation for a single ray. The ray propagates through the model until it reach model edge or a blocking condition is met. The blocking condition depends on the RT calculation phase (see subroutine deposit()).

Parameters
[in]ncID number of the cell originating the ray
[in]alPhi angle of the ray direction
[in]dlTheta angle of the ray direction (ray-direction^Z).
[in]clvlSubdivision level of the cell originating the ray (but not used when the ray doesn't start from the cell nc
[in]src_lumLuminosity associated with the ray
[in]dplane_rayDistance from last intersection plane along axis perpendicular to that plane.
[in]prev_rayDistance crossed by the blocked ray from which the current ray is derived.
[in]isel_rayDirection perpendicular to the last intersection plane (0 = x, 1= y, 2=z).
[in]ray_typeType of ray (e.g. ray_type_high(), ray_type_low(), ray_type_reco(), ray_type_gone(), ray_type_i_obs(), ray_type_i_obs_in() )
[in]ipixHEALPix number of the ray
[in]cc_oldID number of the last cell intersected by the blocked ray (this can be the same or not for the current ray)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::remove_lock_from_cell ( integer  i)

Decreases lock_cell() by a unity.

Here is the caller graph for this function:

subroutine rt_routines::remove_negative ( real(kind=real64), dimension(0:)  array)

Assigns zero to all negative elements of an array. Used in fix_symmetry_part3() to remove negative due to numerical accuracy. Remember that this could give problems in case you want to work with negative intensities.

Here is the caller graph for this function:

subroutine rt_routines::remove_negative_2d ( real(kind=real64), dimension(0:, 0:)  array)

Same as remove_negative() but for 2D arrays.

subroutine rt_routines::rt_loop ( )

This is the main RT loop over the sources of radiation. It is used in all calculation phases where the radiation field is derived.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::rt_loop_2d ( )

Contains the main RT loop over the sources of radiation for the RT 2D mode. The loop is perfomed twice because the first instance processes only the cells with x,y,z > 0 while the second processes those with x =0 OR y=0 OR z=0. Repeating the loops is necessary because only the off-axis cells provide contributions toe the radiation field that can be easily symmetrised.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::rt_loop_iobs ( )

Performs the ray-tracing loop to calculate the output specific intensity arrays i_obs().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::scale_dens_arr ( )

Assigns the values of dens_arr() by scaling the dens() array read from the main grid file.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::scaspe_temp_transfer ( )

Sends chunks of scaspe_temp_send() to the other processes using mpi_scatterv. It stores the received values in scaspe_temp_recv().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine rt_routines::set_en_lim ( )

Sets the en_lim ( $ f_U $) parameter, which checks if a ray carries negligible luminosity and can be blocked without significantly reducing the calculation accuracy.

Here is the caller graph for this function:

subroutine rt_routines::set_i_opacity_arrays ( integer  i0,
integer  i1 
)

Sets the subscript of the first element of the value segment within kext_arr(), gsca_arr() and ksca_norm_arr() used in the RT calculation.

Here is the caller graph for this function:

subroutine rt_routines::set_iq_list ( real(kind=real64), dimension(0:lnum-1)  ray_intensity,
logical  flag_beam 
)

Sets the subscript list iq() of the ray_intensity() elements higher than zero and sets ray_intensity_dep().

Here is the caller graph for this function:

subroutine rt_routines::set_npix_arr ( )

Sets the number of pixels for the scaspe() arrays at each wavelength (npix_arr() and npix_hp_arr()). Having a variable size for the scattering source function allows to save memory especially in the infrared.

Here is the caller graph for this function:

subroutine rt_routines::set_units ( )

Sets units_i_obs(), units_ufield() and cs().

Here is the caller graph for this function:

subroutine rt_routines::set_walls ( )

Sets the values of the wall coordinates in model units.

Here is the caller graph for this function:

subroutine rt_routines::set_wavelength_index ( integer  il,
integer  k 
)

Sets the right wavelength subscript depending whether no_communications() is TRUE or FALSE.

Here is the caller graph for this function:

logical function rt_routines::skip_cell_2d_rt_loop ( integer  i)

Used to skip cell processing within rt_loop_2D. TRUE when the cell i does not have to be processed within 2D loop i2d.

Here is the caller graph for this function:

subroutine rt_routines::sort_en_sca_list ( integer, dimension(0:np_mpi-1)  el_out,
integer, dimension(0:np_mpi-1)  i0_el_out 
)

Sorts the elements of en_sca_list for each thread before passing it to en_sca_list_all(). This avoids to make the sorting later within the master block.

Parameters
el_outNumber of elements to be sent to each MPI process in en_sca_list
i0_el_outAfter the sorting of en_sca_list, this array gives the position of the first element in en_sca_list of the blocks to be sent to each MPI process.

Here is the caller graph for this function:

real(kind = real64) function, dimension(0:lnum-1) rt_routines::src_lum_value ( integer  id)

Assigns the source luminosity within rt_loop() depending on the rt_type() and the cell/point source ID number.

Parameters
[in]idCell/point source ID number

Here is the caller graph for this function:

subroutine rt_routines::store_reshape_arrays ( )

Changes the size of the wavelength dimension of u_fest_arr() to match the dust emission wavelength range. Stores the radiation field energy density in the stellar emission wavelength range into u_final_uv_opt(). New u_final_arr and dens_arr arrays are created to cover the wavelength range of the dust emission. Note that in this subroutine, lnum is still set to lnum_stars.

Here is the caller graph for this function:

subroutine rt_routines::update_i_obs_arr ( )

Adds the i_obs_arr(), calculated durint the current dust heating iteration, to the arrays i_obs_arr_dir() and i_obs_arr_tot().

Here is the caller graph for this function:

subroutine rt_routines::wait_end_thread_loop ( )

Waits until all threads and all MPI processes have completed the source loop and sent out the remaining en_sca elements.

Here is the call graph for this function:

Here is the caller graph for this function:

Variable Documentation

real(kind=real64) rt_routines::accuracy
Parameters
accuracyAccuracy parameter determining value of en_lim() (that is, $ f_U $) through the formula en_lim() = accuracy()/( tot_sources*0.25) with tot_sources = sum(leaf_cells) + tot_p_src() - tot_spare_cells(). Allowed range = [0,1].
type(var_arr_1d), dimension(:), allocatable rt_routines::ads_arr
Parameters
ads_arr(0:dim_npix_unique-1)%a(0:npix_unique-1)Array of cos(theta) used when calculating the values of the Henyey-Greenstein function for a certain ray direction and the angular directions of the scaspe() arrays.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr128
Parameters
ads_arr128(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(2)-1).Precalculated ads factor for nside() = 128.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr16
Parameters
ads_arr16(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(2)-1).Precalculated ads factor for nside() = 16.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr256
Parameters
ads_arr256(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(2)-1).Precalculated ads factor for nside() = 256.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr32
Parameters
ads_arr32(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(2)-1).Precalculated ads factor for nside() = 32.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr4
Parameters
ads_arr4(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(0)-1).Precalculated ads_arr() factor for nside() =4.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr64
Parameters
ads_arr64(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(2)-1).Precalculated ads factor for nside() = 64.
type(var_arr_2d), dimension(:), allocatable rt_routines::ads_arr8
Parameters
ads_arr8(0:dim_npix_unique-1)%a(0:npix_unique-1,0:npix_rays(1)-1).Precalculated ads factor for nside() = 8.
integer(kind=int32) rt_routines::bm_par
Parameters
bm_parMinimum number of ray-intersections per cell during the pre-iteration and direct-light processing for the stellar emission RT. Allowed range = [0,1000].
integer(kind=int32) rt_routines::bm_par_max
Parameters
bm_par_maxMaximum number of ray-intersections per cell. Allowed range = [100*max([bm_par,bm_par_sca]), inf]'
integer(kind=int32) rt_routines::bm_par_sca
Parameters
bm_par_scaMinimum number of ray-intersections per cell during the scattered light processing for the stellar emission RT and all phases for the dust emission RT. Allowed range = [0,1000].
integer, dimension(:,:), allocatable rt_routines::cc_list_all
Parameters
cc_list_all(0np_mpi()-1,0: nproc()-1) ID numbers of the cells whose scaspe_arr() values have to be transmitted between the MPI processes (SHARED)
integer, dimension(:), allocatable rt_routines::cc_list_local
Parameters
cc_list_local(0nproc()-1) ID numbers of the cells whose scaspe_arr() values have to be transmitted to the local MPI processes (SHARED)
integer, dimension(:,:), allocatable rt_routines::ccindd
Parameters
ccindd(3,max_lvl)Tree coordinate position of a cell (usually the intersected cell).
integer, dimension(:,:), allocatable rt_routines::ccindd_nc
Parameters
ccindd_nc(3,max_lvl)Tree coordinate position of the cell originating the ray.
integer rt_routines::chunk_size
Parameters
chunk_sizeChunk size for the OPENMP loops. Set in set_chunk_size().
integer rt_routines::clvl_old
Parameters
clvl_oldSubdivision level of the last intersected cell.
logical rt_routines::cnflag
Parameters
cnflagLogical variable equal to ( tot_lumcell/ tot_rad_en() < conv_en_lim()) witn tot_lumcell=sum( lumcell()). When TRUE, the scattered luminosity left for processing is too small and the scattering iterations are stopped.
logical rt_routines::cnflag_dust
Parameters
cnflag_dustLogical variable used to stop the dust heating iterations. Equal to TRUE when dens_stars_arr() is a small fraction conv_en_lim() of the comulated emission stored in dens_stars_arr_old().
real(kind=real64) rt_routines::conv_en_lim
Parameters
conv_en_limThis is the parameter $ f_L $ used to determine when to stop the scattering iterations. This happens when at the end of an iteration the remaining scattered luminosity to be processed is less than $ f_L $ times the total scattered luminosity (see tot_rad_en()).
type(var_arr_1d), dimension(:), allocatable rt_routines::cos_phi_sca
Parameters
cos_phi_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Array equal to sin( phi_sca()). Used when calculating ads_arr().
type(var_arr_1d), dimension(:), allocatable rt_routines::cos_theta_sca
Parameters
cos_theta_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Array equal to cos( theta_sca()). Used when calculating ads_arr().
integer rt_routines::count_en_sca
Parameters
count_en_scaCounter of elements stored in en_sca_list for the single OPENMP thread (PRIVATE).
integer, dimension(:), allocatable rt_routines::count_en_sca_arr
Parameters
count_en_sca_arrArray of counters of elements stored in en_sca_list for all OPENMP threads (SHARED).
integer rt_routines::count_en_sca_tot
Parameters
count_en_sca_totTotal number of elements to be input in en_sca_list_all() (SHARED).
integer rt_routines::count_transfer
Parameters
count_transferCounter of MPI transfers
integer, parameter rt_routines::count_transfer_start =100
Parameters
count_transfer_startStarting value for count_transfer().
real(kind=real64) rt_routines::cs
Parameters
csLight speed defined in the input.
integer rt_routines::dim_npix_unique
Parameters
dim_npix_uniqueNumber of elements in npix_unique().
real(kind=real64) rt_routines::dist_obs
Parameters
dist_obsObserver distance [pc] used when calculating integrated SED.
logical rt_routines::done
Parameters
doneLogical variable equal to TRUE when the RT algorithm run has been completed.
integer, dimension(:), allocatable rt_routines::el_in_mpi
Parameters
el_in_mpiArray containing the number of en_sca_list element to be received by MPI process (SHARED)
integer, dimension(:,:), allocatable rt_routines::el_out_arr
Parameters
el_out_arr(0nproc()-1, 0: np_mpi() -1) Array containing the number of en_sca_list element to be sent to each MPI process for each OPENMP thread (SHARED).
integer, dimension(:), allocatable rt_routines::el_out_mpi
Parameters
el_out_mpiArray containing the number of en_sca_list element to be sent to each MPI process (SHARED)
real(kind=real64) rt_routines::en_lim
Parameters
en_limThis is the parameter $ f_U $ of the DART-Ray algorithm. It is used when evaluating $ \delta U_\lambda < f_U*U_{\lambda, LL} $
integer rt_routines::en_sca_arrtype
Parameters
en_sca_arrtypeMPI data type ID for en_sca_list_arr() elements SHARED
logical, dimension(:), allocatable rt_routines::en_sca_confirm
Parameters
en_sca_confirmTRUE when the en_sca_list elements from the corresponding MPI process have been received (SHARED)
integer, dimension(:), allocatable rt_routines::en_sca_id_thread
Parameters
en_sca_id_threadContains the ID of the thread that has to process the corresponding element in en_sca_list_received SHARED
type (list_en_sca), dimension(:), allocatable rt_routines::en_sca_list
Parameters
en_sca_listList of scattered luminosities and associated variables that have to be sent to other nodes PRIVATE
type (list_en_sca), dimension(:), allocatable rt_routines::en_sca_list_all
Parameters
en_sca_list_allContains all the en_sca_list() arrays, so the MASTER thread can send this list to the other MPI processes. Note this array is OPEN MP SHARED not PRIVATE. Before the elements are transfered to other arrays, en_sca_list_all() is sorted in order of receiving MPI process.
type (list_en_sca), dimension(:), allocatable rt_routines::en_sca_list_received
Parameters
en_sca_list_receivedContains the scattering parameters received by the local MPI process SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::ffn_arr
Parameters
ffn_arr(0:dim_npix_unique-1)%a(0:npix_unique(i)-1).This is an array used during the calculation of ffn_arr_mw at a single wavelength component. It has to be assigned outside the subroutine to avoid repeated allocation.
type(var_arr_1d), dimension(:), allocatable rt_routines::ffn_arr_mw
Parameters
ffn_arr_mw(0:lnum_node-1)%a(0:npix_arr(i)-1).Scattering phase function values at the angular directions defined for the scaspe() arrays for a given ray direction. A variable length data type is used. The data datype argument contains the values for each line-of-sight direction (npix_arr(i) values for each wavelength i).
real (kind=real64), parameter rt_routines::glepsilon =1.0e-7
Parameters
glepsilonVery small factor used to determine the ray direction (positive or negative). IS THERE A BETTER WAY TO DO THIS ?
logical, dimension(:), allocatable rt_routines::handling_mpi_thread_arr
Parameters
handling_mpi_thread_arrIt is .TRUE. if the corresponding thread has called the handle_mpi_transfer routine SHARED
integer, dimension(:), allocatable rt_routines::i0_count_en_sca
Parameters
i0_count_en_sca(0nproc()* np_mpi()-1) Positions of the first element to be placed within en_sca_list_all by each OPENMP thread for each MPI process (SHARED).
integer, dimension(:), allocatable rt_routines::i0_el_in_mpi
Parameters
i0_el_in_mpiArray containing the first position of the blocks of en_sca_list elements to be received by each MPI process (SHARED).
integer, dimension(:,:), allocatable rt_routines::i0_el_out_arr
Parameters
i0_el_out_arr(0nproc()-1, 0: np_mpi() -1) Array containing the first position of the blocks of en_sca_list elements to be sent to each MPI process for each OPENMP thread (SHARED).
integer, dimension(:), allocatable rt_routines::i0_el_out_mpi
Parameters
i0_el_out_mpiArray containing the first position of the blocks of en_sca_list elements to be sent to each MPI process (SHARED).
integer rt_routines::i2d
Parameters
i2dExternal loop index in rt_loop_2D().
real(kind=real64), dimension(:), allocatable rt_routines::i_obs_temp
Parameters
i_obs_tempTemporary value of specific intensity when calculating i_obs(). It has a different value in each OPENMP thread.
integer, dimension(:), allocatable rt_routines::ik_sca_arr
Parameters
ik_sca_arrFor each wavelength, it gives the right subscript to use when operating with arrays with dimension dim_npix_unique(). For example, it is used to know which element of kp_unique() corresponds to a certain element in kp_sca_arr(). If an element has value -1, it means that scattering is considered isotropic for the corresponding wavelength.
integer, dimension(:), allocatable rt_routines::im_lambda_arr
Parameters
im_lambda_arrArray of the subscripts corresponding to the MPI process which hosts the corresponding wavelengths in the local scaspe arrays (SHARED).
integer, dimension(:), allocatable rt_routines::ind_en_sca_list
Parameters
ind_en_sca_listIndeces of the MPI process to which en_sca_list() elements are transfered PRIVATE
integer, dimension(:), allocatable rt_routines::ind_i_obs
Parameters
ind_i_obsList of indeces of the wavelengths of the i_obs arrays to be printed. By default ind_i_obs contains all the wavelength indeces. If indeces are specified in the input file, only the i_obs for the corresponding wavelengths are printed at the end of the calculation.
integer, dimension(:), allocatable rt_routines::ind_out_maps
Parameters
ind_out_mapsList of indeces of the wavelengths for which the surface brightness maps have to be printed. By default ind_out_maps contains all the wavelength indeces. If indeces are specified in the input file, only the maps for the corresponding wavelengths are printed at the end of the calculation.
integer rt_routines::ipsel_av
Parameters
ipsel_avCounter of rays launched by the same cell. Used in the calculation of the average ray path psel_av_arr().
integer rt_routines::ipsel_av_tot
Parameters
ipsel_av_totCounter of total number of rays
integer, dimension(:), allocatable rt_routines::iq
Parameters
iqSubscript list of the ray_intensity element higher than zero. This array varies its size during the calculation. PRIVATE
logical, dimension(:), allocatable rt_routines::iq_a
Parameters
iq_aArray whose elements are equal to .TRUE. if the corresponding ray_intensity element is higher than zero. This array does not vary its size during the calculation. PRIVATE
integer, dimension(:), allocatable rt_routines::iq_maps_id
Parameters
iq_maps_idArray of the subscripts corresponding to the wavelength index of the locally stored i_obs_arr() arrays for which the surface brightness maps have to be calculated (SHARED).
integer, dimension(:), allocatable rt_routines::iq_sca_id
Parameters
iq_sca_idArray of the subscripts corresponding to the wavelengths stored in the local scaspe arrays (SHARED).
logical, dimension(:), allocatable rt_routines::iq_sca_node
Parameters
iq_sca_nodeArrays whose elements are .TRUE. if corresponding wavelength in lambda_arr() is considered in the local scaspe and i_obs arrays. SHARED
integer rt_routines::iscaspe_big
Parameters
iscaspe_bigCounter of used scaspe_temp_arr_big() elements
integer rt_routines::iterations
Parameters
iterationsCounter of scattering iterations.
integer rt_routines::iterations_dustem
Parameters
iterations_dustemCounter of dust self-heating iterations.
type(var_arr_1d), dimension(:), allocatable rt_routines::ix
Parameters
ix(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after X symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::ixy
Parameters
ixy(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after XY symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::ixyz
Parameters
ixyz(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after XYZ symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::ixz
Parameters
ixz(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after XZ symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::iy
Parameters
iy(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after Y symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::iyz
Parameters
iyz(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after YZ symmetry transformation. SHARED
type(var_arr_1d), dimension(:), allocatable rt_routines::iz
Parameters
iz(0:dim_npix_unique-1)%a(0:npix_hp_unique(i)-1)Scaspe() array indices after Z symmetry transformation. SHARED
integer, dimension(:), allocatable rt_routines::kp_sca_arr
Parameters
kp_sca_arrValues of kp parameter used to determine the values of npix_arr() and npix_hp_arr() at each wavelength.
integer(kind=int32) rt_routines::kp_sca_max
Parameters
kp_sca_maxHEALPix parameter used to determine the maximum allowed number of "spherical pixels" in the scaspe() arrays: npix_max = 12*2^(2*kp_sca_max). The allowed range is [0,4]. Values higher than 4 are not allowed because the scaspe arrays would become too big. However, this can be changed in check_input().
integer, dimension(:), allocatable rt_routines::kp_unique
Parameters
kp_uniqueList of HEALPix kp values corresponding to the values in npix_unique (remember that npix contains the number of HEALPix pixel plus the external observer lines of sight).
integer rt_routines::lnum_a
Parameters
lnum_aNumber of ray_intensity elements higher than zero
integer rt_routines::lnum_a_old
Parameters
lnum_a_oldNumber of ray_intensity elements higher than zero in the previous ray-cell intersection step
integer rt_routines::lnum_maps
Parameters
lnum_mapsNumber of wavelengths for which the surface brightness maps have to be calculated.
integer rt_routines::lnum_node
Parameters
lnum_nodeNumber of wavelengths considered in the scaspe and i_obs arrays stored within the node memory.
integer, dimension(:), allocatable rt_routines::lnum_node_arr
Parameters
lnum_node_arrArray of number of wavelengths considered in the scaspe and i_obs arrays in each node memory. SHARED
integer rt_routines::lnum_node_maps
Parameters
lnum_node_mapsNumber of wavelengths for which the surface brightness maps have to be calculated and that are stored in the local i_obs arrays.
real(kind=real64), dimension(:), allocatable rt_routines::lum_lost
Parameters
lum_lostTotal luminosity not processed in the entire RT run. SHARED
real(kind=real64), dimension(:), allocatable rt_routines::lum_lost_prev
Parameters
lum_lost_prevTotal luminosity not processed during the direct-light processing or the previous scattering iteration. SHARED
real(kind=real64), dimension(:), allocatable rt_routines::lum_lost_temp
Parameters
lum_lost_tempTemporary value of the luminosity not processed within a HEALPix section (currently corresponding to nside = 4). These values accumulate on lum_lost().
integer, parameter rt_routines::max_n_input = 1000
Parameters
max_n_inputMaximum number of input elements for allocatable input arrays as ind_i_obs() and ind_out_maps().
integer rt_routines::max_npix
Parameters
max_npixMaximum value of npix_arr() within the wavelength range used in the current RT algorithm.
logical, dimension(:), allocatable rt_routines::mpi_proc_completed
Parameters
mpi_proc_completedIt is .TRUE. if the corresponding MPI process does not have any other en_sca element to send out. If all its threads are in the waiting part at the end of the source loop in rt_loop or rt_loop_2d, it is used to exit the loop when all elements are .TRUE. SHARED
integer*1, parameter rt_routines::mpi_transfer_en_sca = 1
Parameters
MPI_TRANSFER_EN_SCAID used in handle_mpi_transfers() to distinguish call requiring process
integer*1, parameter rt_routines::mpi_transfer_scaspe = 0
Parameters
MPI_TRANSFER_SCASPEID used in handle_mpi_transfers() to distinguish call requiring transfer of scaspe arrays
integer rt_routines::n_el_in_tot
Parameters
n_el_in_totTotal number of en_sca_list_received elements
logical rt_routines::no_communications
Parameters
no_communicationsTRUE if no communication mode has to be used. This means that the scaspe and i_obs arrays are not distributed among the different MPI processes. Instead, copies of the entire arrays are created for each MPI process. This mode requires much more memory for each node and thus it is suggested only when very few wavelengths are used or the 3D grid resolution is very small. The advantage of using this method is that no communication is performed during the RT calculation. At the end of each stage, the arrays are reduced as for the u_final() arrays. This mode requires sequential_scattering() = .TRUE. because otherwise there is a complex race condition to avoid in assign_scaspe_temp_arr().
integer, dimension(:), allocatable rt_routines::npix_arr
Parameters
npixNumber of pixels in the scaspe() arrays, including user - defined directions (see file_dir_out())
npix_arrNumber of pixels in the scaspe() arrays at each wavelength, including user - defined directions (see file_dir_out())
integer rt_routines::npix_hp
Parameters
npix_hpNumber of pixels in the scaspe() arrays, NOT including user_defined directions
integer, dimension(:), allocatable rt_routines::npix_hp_arr
Parameters
npix_hp_arrNumber of pixels in the scaspe() arrays at each wavelength, NOT including user_defined directions
integer, dimension(:), allocatable rt_routines::npix_hp_unique
Parameters
npix_hp_uniqueList of the single values present in npix_hp_arr() without repetitions.
integer, parameter rt_routines::npix_main =12
Parameters
npix_mainNumber of spherical pixels in the HEALPix sphere for nside()=1.
integer, dimension(:), allocatable rt_routines::npix_unique
Parameters
npix_uniqueList of the single values present in npix_arr() without repetitions.
integer rt_routines::nproc
Parameters
nprocNumber of processors to be used for OPENMP loops.
integer rt_routines::nside
Parameters
nsideHEALPix nside parameter.
integer, parameter rt_routines::nside_min = 4
Parameters
nside_minHEALPix minimum nside parameter used to select the ray angular density in main_dir_loop().
integer(kind=int32) rt_routines::nside_sca
Parameters
nside_scaHEALPix nside parameter. Used when calling some HEALPix subroutines.
integer rt_routines::num_processed_cells
Parameters
num_processed_cellsCounter of processed cells
integer rt_routines::num_scaspe_pass
Parameters
num_scaspe_passNumber of cells whose scaspe values are transfered to other nodes.
logical rt_routines::only_direct_rt
Parameters
only_direct_rtTRUE if only direct light has to be processed.
real(kind=real64), dimension(3) rt_routines::pabs_arr
Parameters
pabs_arr(3)Array storing some constant factors when calculating the ray-cell intersection in calc_psel(). These factors are the same as in the previous intersected cell if the subdivision level, stored in clvl_old(), does not change.
integer(kind=int32) rt_routines::pabs_max
Parameters
pabs_maxMaximum value allowed for pabs variable in calc_psel(). Set to 10* modelsize()
character(len=lcar) rt_routines::param_to_project
Parameters
param_to_projectParameter to be projected in the 'projection' rt_algorithm(). Choices are: 'stellar_emission' and 'optical_depth'.
type(var_arr_1d), dimension(:), allocatable rt_routines::phi_sca
Parameters
phi_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Arrays of phi angles (n^xy) corresponding to the elements of the first index of the scaspe() arrays.
logical rt_routines::print_psel_av
Parameters
print_psel_avEqual to TRUE if psel_av_arr() array has to be in the output.
logical rt_routines::print_scaspe_tot
Parameters
print_scaspe_totLogical equal to TRUE if scaspe_tot() has to be printed at the end of the scattering iterations. The scaspe_tot() array can then be used in RT algorithm "i_obs" (see rt_algorithm()).
real(kind=real64) rt_routines::psel_av
Parameters
psel_avAverage ray path for the rays cast from a certain cell.
real(kind=real64) rt_routines::rad_lim
Parameters
rad_limMaximum distance relative to modelsize() that can be crossed by the rays in the precalculation phase. Allowed range = [0,2].
integer, parameter rt_routines::ras_first_launch = 3
Parameters
ras_first_launchray_status ID value for rays that have just been cast from the corresponding radiation source.
integer, parameter rt_routines::ras_go_high = 0
Parameters
ras_go_highray_status ID value for rays that have to be split.
integer, parameter rt_routines::ras_go_low = 2
Parameters
ras_go_lowray_status ID value for rays that can be potentially merged.
integer, parameter rt_routines::ras_gone = 1
Parameters
ras_goneray_status ID value for rays that can be eliminated.
integer, parameter rt_routines::ras_re_launched = 4
Parameters
ras_re_launchedray_status ID value for rays that are child/parent rays and that start their propagation where the corresponding parent/child rays have been blocked
real(kind=real64), dimension(:), allocatable rt_routines::ray_intensity_dep
Parameters
ray_intensity_depArray containing the elements of ray_intensity higher than zero. It is input to the subroutine deposit (you cannot input ray_intensity(iq) into a subroutine. The array elements are not updated.) PRIVATE
character(len = 20) rt_routines::rt_algorithm
Parameters
rt_algorithmType of RT algorithm to be used. Choices are: 'main': this is the standard RT algorithm. Firstly, it proceeds with the calculation of the radiation fields, scattering source function, outgoing radiation specific intensity, surface brightness maps for the stellar emission. Then, it performs the same for the dust emission (unless no_dust_rt() = .TRUE.); '2D': same as the 'main' RT algorithm but for axisymmetric models. In this case, the code checks that the input grid is indeed axisymmetric and then performs the RT calculations taking into account the symmetries. This mode performs RT calculations up to a factor 8 faster then the main mode; 'dust' and 'dust_2D': these RT algorithms perform the dust emission RT calculations. To use these modes, it is required that the stellar emission RT has already been performed and the output radiation fields are printed on disk. The two modes are for the normal and 2D algorithm respectively; 'sed' and 'sed_dust': these modes can be used to calculate the integrated SEDs and the maps using the outgoing specific intensity files ( file_i_obs() and file_i_obs_in()) printed in a previous RT calculation. The two modes are for the stellar emission and the dust emission respectively; 'i_obs' and 'i_obs_dust': these modes can be used to calculate the outgoing specific intensity and the corresponding surface brightness maps for arbitrary observer positions without repeating a complete RT calculation. These modes can be used only if the scattering source function has be printed on disk (see print_scaspe_tot()). The two modes are for the stellar emission and dust emission respectively; 'projection': this mode can be used to make maps of the stellar emission without any dust or maps of the dust optical depth. Useful to check that the input grid has been correctly calculated before starting a complete RT calculation.
integer rt_routines::rt_algorithm_id
Parameters
rt_algorithm_IDInteger variable associated with the rt_algorithm() choices. See check_input() for definitions.
integer rt_routines::rt_type
Parameters
rt_typeInteger variable used to select the different types of RT algorithm steps (e.g. precalculation or direct light processing). See select_rt_type().
integer, parameter rt_routines::rta_2d = 2
Parameters
rta_2Drt_algorithm_ID() ID value for the 2D RT calculation (still 3D but it exploits the 2D symmetry of the input model)
integer, parameter rt_routines::rta_dust = 4
Parameters
rta_dustrt_algorithm_ID() ID value for the dust RT calculation.
integer, parameter rt_routines::rta_dust2d = 5
Parameters
rta_dust_2Drt_algorithm_ID() ID value for the dust RT calculation in 2D mode.
integer, parameter rt_routines::rta_i_obs = 1
Parameters
rta_i_obsrt_algorithm_ID() ID value for the outgoing radiation specific intensity RT calculation for the stellar emission.
integer, parameter rt_routines::rta_i_obs_dust = 6
Parameters
rta_i_obs_dustrt_algorithm_ID() ID value for the outgoing radiation specific intensity RT calculation for the dust emission.
integer, parameter rt_routines::rta_main = 0
Parameters
rta_mainrt_algorithm_ID() ID value for the standard RT calculation.
integer, parameter rt_routines::rta_projection = 8
Parameters
rta_projectionrt_algorithm_ID() ID value for projection of physical quantities as optical depth, intrinsic stellar emission.
integer, parameter rt_routines::rta_sed = 3
Parameters
rta_sedrt_algorithm_ID() ID value for the SED calculation.
integer, parameter rt_routines::rta_sed_dust = 7
Parameters
rta_sed_dustrt_algorithm_ID() ID value for the SED calculation for the dust emission.
integer, parameter rt_routines::rtt_dir_cell = 5
Parameters
rtt_dir_cellrt_type() ID for the direct light calculation for emitting cells.
integer, parameter rt_routines::rtt_dir_src = 7
Parameters
rtt_dir_srcrt_type() ID for the direct light calculation for point sources.
integer, parameter rt_routines::rtt_grid_init_dust = 17
Parameters
rtt_grid_init_dustrt_type() ID for dust emission RT grid initilialization.
integer, parameter rt_routines::rtt_grid_init_projection = 18
Parameters
rtt_grid_init_projectionID for projection algorithm grid initilialization.
integer, parameter rt_routines::rtt_grid_init_stars = 16
Parameters
rtt_grid_init_starsrt_type() ID for stellar emission RT grid initilialization.
integer, parameter rt_routines::rtt_i_obs = 11
Parameters
rtt_i_obsrt_type() ID for the outgoing radiation specific intensity calculation for scattered radiation sources.
integer, parameter rt_routines::rtt_i_obs_dir_cell = 6
Parameters
rtt_i_obs_dir_cellrt_type() ID for the outgoing radiation specific intensity calculation for emitting cells.
integer, parameter rt_routines::rtt_i_obs_dir_src = 8
Parameters
rtt_i_obs_dir_srcrt_type() ID for the outgoing radiation specific intensity calculation for point sources.
integer, parameter rt_routines::rtt_output_part1 = 3
Parameters
rtt_output_part1rt_type() ID for the output of the results of the precalculation.
integer, parameter rt_routines::rtt_output_part2 = 9
Parameters
rtt_output_part2rt_type() ID for the output of the direct light processing.
integer, parameter rt_routines::rtt_precalc_cell = 1
Parameters
rtt_precalc_cellrt_type() ID for the precalculation of the radiation field due to the emitting cells.
integer, parameter rt_routines::rtt_precalc_src = 2
Parameters
rtt_precalc_srcrt_type() ID for the precalculation of the radiation field due to the point sources.
integer, parameter rt_routines::rtt_prep_part2 = 4
Parameters
rtt_prep_part2rt_type() ID for the preparation of the direct light processing.
integer, parameter rt_routines::rtt_read_i_obs = 13
Parameters
rtt_read_i_obsrt_type() ID for the reading of i_obs() array (direct light and scattered light included) from file.
integer, parameter rt_routines::rtt_read_i_obs_part2 = 12
Parameters
rtt_read_i_obs_part2rt_type() ID for the reading of i_obs() array (direct light only) from file.
integer, parameter rt_routines::rtt_read_scaspe_tot = 15
Parameters
rtt_read_scaspe_totrt_type() ID for the reading of the scaspe_tot() array from file.
integer, parameter rt_routines::rtt_read_ufield = 14
Parameters
rtt_read_ufieldrt_type() ID for the reading of u_field() array from file.
integer, parameter rt_routines::rtt_scatt = 10
Parameters
rtt_scattrt_type() ID for the scattered radiation calculation.
integer, parameter rt_routines::rtt_start = 0
Parameters
rtt_startrt_type() ID for the beginning of the RT run.
type(var_arr_1d), dimension(:), allocatable rt_routines::scaspe_temp_arr
Parameters
scaspe_temp_arr(0:lnum-1)%a(0npix_hp_arr()+ tot_ndir_scaspe()-1) Extracted part of the scaspe_arr() array used by a single thread. (PRIVATE)
type(var_arr_2d), dimension(:), allocatable rt_routines::scaspe_temp_arr_big
Parameters
scaspe_temp_arr_big(0:lnum-1)%a(0npix_hp_arr()+ tot_ndir_scaspe()-1,0:num_scaspe_pass-1) Storage of the extracted parts of the scaspe_arr() array used by a single thread. This is used so the scaspe values are transfered less often (PRIVATE)
type(var_arr_2d), dimension(:), allocatable rt_routines::scaspe_temp_recv
Parameters
scaspe_temp_recv(0:lnum-1)%a(0npix_hp()+ tot_ndir_scaspe()-1,0: nproc()*num_scaspe_pass()-1) Scaspe_arr() elements received by the local MPI process (SHARED)
real(kind=real64), dimension(:), allocatable rt_routines::scaspe_temp_send
Parameters
scaspe_temp_send(0sum(npix_arr(i))* (nproc()* np_mpi()*num_scaspe_pass)-1) with the sum over npix_arr only for i corresponding to the wavelengths locally stored in the scaspe_arr() array. Scaspe_arr() elements sent by the local MPI process (SHARED).
real(kind=real64), dimension(:,:), allocatable rt_routines::sed_arr
Parameters
sed_arr(0lnum_tot()-1, 0: tot_ndir() -1) Contains total emission SED including all components and for all directions (stellar or dust emission depending on the algorithm).
real(kind=real64), dimension(:,:), allocatable rt_routines::sed_arr_dir
Parameters
sed_arr_dir(0lnum_tot()-1, 0: tot_ndir() -1) Contains emission SED of the direct light (stellar or dust emission depending on the algorithm).
logical rt_routines::sequential_scattering
Parameters
sequential_scatteringTrue if, during each scattering iterations thescattered luminosity is propagated for the same scattering order. That is, the contribution to the scattered radiation, calculated after each ray-cell intersection, are stored in scaspe_arr but the assignment of scaspe_temp is done using scaspe_prev_arr which contains the scaspe_arr calculated in the previous iteration (before the current iteration starts, scaspe_arr is initialised).
type(var_arr_1d), dimension(:), allocatable rt_routines::sin_phi_sca
Parameters
sin_phi_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Array equal to sin( phi_sca()). Used when calculating ads_arr().
type(var_arr_1d), dimension(:), allocatable rt_routines::sin_theta_sca
Parameters
sin_theta_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Array equal to sin(theta_sca()). Used when calculating ads_arr().
integer rt_routines::size_en_sca_list
Parameters
size_en_sca_listSize of en_sca_list().
real(kind=real64), dimension(3) rt_routines::src_ccoord
Parameters
src_ccoord(3)Radiation source 3D position (it can be either a cell centre or a point source position).
real(kind=real64) rt_routines::tau_cell_max
Parameters
tau_cell_maxMaximum optical depth allowed for leaf cells in the grid creation program. Allowed range = [0,0.5].
type (list_en_sca), dimension(:), allocatable rt_routines::temp_en_sca_list
Parameters
temp_en_sca_listCopy of en_sca_list() used when sorting array PRIVA
logical rt_routines::test_run
Parameters
test_runIf set TRUE it will not run the rt_loop routines. Useful to check that the all the subroutines are running fine and the output arrays are correctly printed. However, the results will not be correct and therefore should be deleted afterwards.
type(var_arr_1d), dimension(:), allocatable rt_routines::theta_sca
Parameters
theta_sca(0:dim_npix_unique-1)%a(0:npix_unique(i)-1)Array of theta angles (n^z) corresponding to the elements of the first index of the scaspe() arrays.
integer rt_routines::threadid
Parameters
ThreadIDOPENMP thread ID number
real(kind=real64), parameter rt_routines::tol_p = 1E-6

param tol_p_src Tolerance parameter to find host cell of a point source

integer, dimension(:), allocatable rt_routines::tot_npix_arr_local
Parameters
tot_npix_arr_localTotal number of elements in npix_arr for the wavelength locally stored in the scaspe_arr() array. Used to calculate how many elements have to be stored in scaspe_temp_send().
real(kind=real64), dimension(:), allocatable rt_routines::tot_rad_en
Parameters
tot_rad_enTotal luminosity for all wavelengths. During precalculation and direct-light processing is the total luminosity of the stars. During the scattering iteration is the total scattered luminosity.
real(kind=real64), dimension(:), allocatable rt_routines::tot_rad_en_or
Parameters
tot_rad_en_orTotal stellar luminosity for all wavelengths. Used to determine lost luminosity fraction in print_lum_lost(). If this fraction is higher than 1%, a warning is output on the terminal.
integer*1, dimension(:,:), allocatable rt_routines::transfer_type_all
Parameters
transfer_type_all(0nproc()-1, 0: np_mpi()-1) Array containing the transfer_type value of each OPENMP thread and MPI process SHARED
integer*1, dimension(:), allocatable rt_routines::transfer_type_local
Parameters
transfer_type_local(0nproc()-1) Array containing the transfer_type values for the local OPENMP threads. SHARED
integer rt_routines::transfer_type_tot_all
Parameters
transfer_type_tot_allSum of transfer_type_tot_local()
integer*1 rt_routines::transfer_type_tot_local
Parameters
transfer_type_tot_localSum of transfer_type_local()
real(kind=real64) rt_routines::vec_mod
Parameters
vec_modDistance between radiation source and the internal observer.
logical, dimension(:), allocatable rt_routines::wait_thread_arr
Parameters
wait_thread_arrIt is .TRUE. if the corresponding thread has entered the waiting region at the end of the source loop in rt_loop or rt_loop_2d. SHARED
real(kind=real64), dimension(2) rt_routines::x_wall_coord
Parameters
x_wall_coordPositions of the X-perpendicular walls in relative coordinates (0 = -modelsize/2., 1 = modelsize/2).
logical rt_routines::x_wall_on
Parameters
x_wall_onTRUE if wall perpendicular to X direction is set
real(kind=real64), dimension(2) rt_routines::y_wall_coord
Parameters
y_wall_coordPositions of the Y-perpendicular walls in relative coordinates (0 = -modelsize/2., 1 = modelsize/2).
logical rt_routines::y_wall_on
Parameters
y_wall_onTRUE if wall perpendicular to Y direction is set
real(kind=real64), dimension(2) rt_routines::z_wall_coord
Parameters
z_wall_coordPositions of the Z-perpendicular walls in relative coordinates (0 = -modelsize/2., 1 = modelsize/2).
logical rt_routines::z_wall_on
Parameters
z_wall_onTRUE if wall perpendicular to Z direction is set