classes_neighbour_cells.f90 Source File

This File Depends On

sourcefile~~classes_neighbour_cells.f90~~EfferentGraph sourcefile~classes_neighbour_cells.f90 classes_neighbour_cells.f90 sourcefile~classes_pair_potential.f90 classes_pair_potential.f90 sourcefile~classes_pair_potential.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~classes_visitable_walls.f90 classes_visitable_walls.f90 sourcefile~classes_pair_potential.f90->sourcefile~classes_visitable_walls.f90 sourcefile~classes_parallelepiped_domain.f90 classes_parallelepiped_domain.f90 sourcefile~classes_parallelepiped_domain.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~procedures_checks.f90 procedures_checks.f90 sourcefile~procedures_checks.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~procedures_checks.f90->sourcefile~classes_pair_potential.f90 sourcefile~procedures_checks.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~classes_dipolar_neighbourhood.f90 classes_dipolar_neighbourhood.f90 sourcefile~procedures_checks.f90->sourcefile~classes_dipolar_neighbourhood.f90 sourcefile~classes_potential_expression.f90 classes_potential_expression.f90 sourcefile~procedures_checks.f90->sourcefile~classes_potential_expression.f90 sourcefile~procedures_checks.f90->sourcefile~classes_visitable_walls.f90 sourcefile~classes_periodic_box.f90 classes_periodic_box.f90 sourcefile~procedures_checks.f90->sourcefile~classes_periodic_box.f90 sourcefile~classes_floor_penetration.f90 classes_floor_penetration.f90 sourcefile~procedures_checks.f90->sourcefile~classes_floor_penetration.f90 sourcefile~classes_min_distance.f90 classes_min_distance.f90 sourcefile~procedures_checks.f90->sourcefile~classes_min_distance.f90 sourcefile~classes_dirac_distribution_plus.f90 classes_dirac_distribution_plus.f90 sourcefile~procedures_checks.f90->sourcefile~classes_dirac_distribution_plus.f90 sourcefile~data_cells.f90 data_cells.f90 sourcefile~data_cells.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~classes_hard_contact.f90 classes_hard_contact.f90 sourcefile~classes_hard_contact.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~data_constants.f90 data_constants.f90 sourcefile~data_constants.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~data_constants.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~data_constants.f90->sourcefile~procedures_checks.f90 sourcefile~data_constants.f90->sourcefile~data_cells.f90 sourcefile~data_constants.f90->sourcefile~classes_visitable_walls.f90 sourcefile~procedures_parallelepiped_domain_micro.f90 procedures_parallelepiped_domain_micro.f90 sourcefile~data_constants.f90->sourcefile~procedures_parallelepiped_domain_micro.f90 sourcefile~data_constants.f90->sourcefile~classes_periodic_box.f90 sourcefile~data_constants.f90->sourcefile~classes_floor_penetration.f90 sourcefile~data_constants.f90->sourcefile~classes_dirac_distribution_plus.f90 sourcefile~procedures_dipolar_interactions_micro.f90 procedures_dipolar_interactions_micro.f90 sourcefile~data_constants.f90->sourcefile~procedures_dipolar_interactions_micro.f90 sourcefile~classes_dipolar_neighbourhood.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~procedures_errors.f90 procedures_errors.f90 sourcefile~procedures_errors.f90->sourcefile~classes_neighbour_cells.f90 sourcefile~procedures_errors.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~procedures_errors.f90->sourcefile~procedures_checks.f90 sourcefile~procedures_errors.f90->sourcefile~classes_visitable_walls.f90 sourcefile~procedures_errors.f90->sourcefile~classes_periodic_box.f90 sourcefile~procedures_errors.f90->sourcefile~classes_floor_penetration.f90 sourcefile~classes_potential_expression.f90->sourcefile~classes_pair_potential.f90 sourcefile~types_potential_domain.f90 types_potential_domain.f90 sourcefile~types_potential_domain.f90->sourcefile~classes_pair_potential.f90 sourcefile~types_potential_domain.f90->sourcefile~procedures_checks.f90 sourcefile~types_potential_domain_selector.f90 types_potential_domain_selector.f90 sourcefile~types_potential_domain_selector.f90->sourcefile~classes_pair_potential.f90 sourcefile~types_potential_domain_selector.f90->sourcefile~procedures_checks.f90 sourcefile~classes_visitable_walls.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~procedures_parallelepiped_domain_micro.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~classes_periodic_box.f90->sourcefile~classes_parallelepiped_domain.f90 sourcefile~classes_periodic_box.f90->sourcefile~classes_visitable_walls.f90 sourcefile~classes_floor_penetration.f90->sourcefile~classes_visitable_walls.f90 sourcefile~classes_min_distance.f90->sourcefile~classes_visitable_walls.f90 sourcefile~procedures_centered_block_micro.f90 procedures_centered_block_micro.f90 sourcefile~procedures_centered_block_micro.f90->sourcefile~classes_floor_penetration.f90 sourcefile~classes_number_to_string.f90 classes_number_to_string.f90 sourcefile~classes_number_to_string.f90->sourcefile~procedures_checks.f90 sourcefile~data_strings.f90 data_strings.f90 sourcefile~data_strings.f90->sourcefile~classes_number_to_string.f90 sourcefile~classes_dirac_distribution_plus.f90->sourcefile~classes_hard_contact.f90 sourcefile~procedures_dipolar_interactions_micro.f90->sourcefile~classes_dipolar_neighbourhood.f90
Help

Files Dependent On This One

sourcefile~~classes_neighbour_cells.f90~~AfferentGraph sourcefile~classes_neighbour_cells.f90 classes_neighbour_cells.f90 sourcefile~classes_visitable_cells.f90 classes_visitable_cells.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~classes_visitable_cells.f90 sourcefile~procedures_short_interactions_resetter.f90 procedures_short_interactions_resetter.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~procedures_short_interactions_resetter.f90 sourcefile~procedures_neighbour_cells_factory.f90 procedures_neighbour_cells_factory.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~procedures_neighbour_cells_factory.f90 sourcefile~classes_visitable_cells_memento.f90 classes_visitable_cells_memento.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~classes_visitable_cells_memento.f90 sourcefile~procedures_visitable_cells_factory.f90 procedures_visitable_cells_factory.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~procedures_visitable_cells_factory.f90 sourcefile~types_cells_wrapper.f90 types_cells_wrapper.f90 sourcefile~classes_neighbour_cells.f90->sourcefile~types_cells_wrapper.f90 sourcefile~classes_visitable_cells.f90->sourcefile~procedures_short_interactions_resetter.f90 sourcefile~classes_visitable_cells.f90->sourcefile~classes_visitable_cells_memento.f90 sourcefile~classes_visitable_cells.f90->sourcefile~procedures_visitable_cells_factory.f90 sourcefile~classes_visitable_cells.f90->sourcefile~types_cells_wrapper.f90 sourcefile~classes_dipolar_neighbourhoods_visitor.f90 classes_dipolar_neighbourhoods_visitor.f90 sourcefile~classes_visitable_cells.f90->sourcefile~classes_dipolar_neighbourhoods_visitor.f90 sourcefile~procedures_short_interactions_visitor.f90 procedures_short_interactions_visitor.f90 sourcefile~classes_visitable_cells.f90->sourcefile~procedures_short_interactions_visitor.f90 sourcefile~classes_maximum_box_compression_explorer.f90 classes_maximum_box_compression_explorer.f90 sourcefile~classes_visitable_cells.f90->sourcefile~classes_maximum_box_compression_explorer.f90 sourcefile~classes_boxes_volume_exchange.f90 classes_boxes_volume_exchange.f90 sourcefile~procedures_short_interactions_resetter.f90->sourcefile~classes_boxes_volume_exchange.f90 sourcefile~classes_volume_change_method.f90 classes_volume_change_method.f90 sourcefile~procedures_short_interactions_resetter.f90->sourcefile~classes_volume_change_method.f90 sourcefile~classes_box_volume_change.f90 classes_box_volume_change.f90 sourcefile~procedures_short_interactions_resetter.f90->sourcefile~classes_box_volume_change.f90 sourcefile~procedures_plmc_resetter.f90 procedures_plmc_resetter.f90 sourcefile~procedures_short_interactions_resetter.f90->sourcefile~procedures_plmc_resetter.f90 sourcefile~procedures_visitable_cells_memento_factory.f90 procedures_visitable_cells_memento_factory.f90 sourcefile~classes_visitable_cells_memento.f90->sourcefile~procedures_visitable_cells_memento_factory.f90 sourcefile~procedures_cells_memento.f90 procedures_cells_memento.f90 sourcefile~classes_visitable_cells_memento.f90->sourcefile~procedures_cells_memento.f90 sourcefile~types_short_interactions_wrapper.f90 types_short_interactions_wrapper.f90 sourcefile~classes_visitable_cells_memento.f90->sourcefile~types_short_interactions_wrapper.f90 sourcefile~procedures_visitable_cells_factory.f90->sourcefile~classes_visitable_cells_memento.f90 sourcefile~types_cells_wrapper.f90->sourcefile~classes_boxes_volume_exchange.f90 sourcefile~types_cells_wrapper.f90->sourcefile~classes_volume_change_method.f90 sourcefile~types_cells_wrapper.f90->sourcefile~classes_box_volume_change.f90 sourcefile~types_cells_wrapper.f90->sourcefile~procedures_cells_memento.f90 sourcefile~types_cells_wrapper.f90->sourcefile~types_short_interactions_wrapper.f90 sourcefile~procedures_exchange_visitors.f90 procedures_exchange_visitors.f90 sourcefile~types_cells_wrapper.f90->sourcefile~procedures_exchange_visitors.f90 sourcefile~procedures_transmutation_updaters.f90 procedures_transmutation_updaters.f90 sourcefile~types_cells_wrapper.f90->sourcefile~procedures_transmutation_updaters.f90 sourcefile~procedures_exchange_updaters.f90 procedures_exchange_updaters.f90 sourcefile~types_cells_wrapper.f90->sourcefile~procedures_exchange_updaters.f90 sourcefile~procedures_transmutation_visitors.f90 procedures_transmutation_visitors.f90 sourcefile~types_cells_wrapper.f90->sourcefile~procedures_transmutation_visitors.f90 sourcefile~types_markov_chain_explorer_wrapper.f90 types_markov_chain_explorer_wrapper.f90 sourcefile~classes_dipolar_neighbourhoods_visitor.f90->sourcefile~types_markov_chain_explorer_wrapper.f90 sourcefile~procedures_dipolar_neighbourhoods_visitors_factory.f90 procedures_dipolar_neighbourhoods_visitors_factory.f90 sourcefile~classes_dipolar_neighbourhoods_visitor.f90->sourcefile~procedures_dipolar_neighbourhoods_visitors_factory.f90 sourcefile~procedures_exploration_inquirers.f90 procedures_exploration_inquirers.f90 sourcefile~classes_dipolar_neighbourhoods_visitor.f90->sourcefile~procedures_exploration_inquirers.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~classes_dipolar_neighbourhoods_visitor.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~classes_maximum_box_compression_explorer.f90 sourcefile~procedures_plmc_visitor.f90 procedures_plmc_visitor.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~procedures_plmc_visitor.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~classes_boxes_volume_exchange.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~classes_volume_change_method.f90 sourcefile~procedures_short_interactions_visitor.f90->sourcefile~classes_box_volume_change.f90 sourcefile~classes_maximum_box_compression_explorer.f90->sourcefile~types_markov_chain_explorer_wrapper.f90 sourcefile~procedures_maximum_boxes_compression_explorer_factory.f90 procedures_maximum_boxes_compression_explorer_factory.f90 sourcefile~classes_maximum_box_compression_explorer.f90->sourcefile~procedures_maximum_boxes_compression_explorer_factory.f90 sourcefile~procedures_markov_chain_explorer_factory.f90 procedures_markov_chain_explorer_factory.f90 sourcefile~types_markov_chain_explorer_wrapper.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_exploring_writers_factory.f90 procedures_exploring_writers_factory.f90 sourcefile~types_markov_chain_explorer_wrapper.f90->sourcefile~procedures_exploring_writers_factory.f90 sourcefile~plmc_explore.f90 plmc_explore.f90 sourcefile~types_markov_chain_explorer_wrapper.f90->sourcefile~plmc_explore.f90 sourcefile~procedures_dipolar_neighbourhoods_visitors_factory.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_exploration_inquirers.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_exploration_inquirers.f90->sourcefile~procedures_exploring_writers_factory.f90 sourcefile~procedures_exploration_inquirers.f90->sourcefile~plmc_explore.f90 sourcefile~procedures_dipolar_interactions_facades_factory.f90 procedures_dipolar_interactions_facades_factory.f90 sourcefile~procedures_exploration_inquirers.f90->sourcefile~procedures_dipolar_interactions_facades_factory.f90 sourcefile~procedures_short_interactions_factory.f90 procedures_short_interactions_factory.f90 sourcefile~procedures_exploration_inquirers.f90->sourcefile~procedures_short_interactions_factory.f90 sourcefile~procedures_markov_chain_explorer_factory.f90->sourcefile~plmc_explore.f90 sourcefile~procedures_physical_model_factory.f90 procedures_physical_model_factory.f90 sourcefile~procedures_dipolar_interactions_facades_factory.f90->sourcefile~procedures_physical_model_factory.f90 sourcefile~procedures_short_interactions_factory.f90->sourcefile~procedures_physical_model_factory.f90 sourcefile~procedures_plmc_visitor.f90->sourcefile~plmc_explore.f90 sourcefile~plmc_generate.f90 plmc_generate.f90 sourcefile~procedures_plmc_visitor.f90->sourcefile~plmc_generate.f90 sourcefile~procedures_boxes_volume_exchange_factory.f90 procedures_boxes_volume_exchange_factory.f90 sourcefile~classes_boxes_volume_exchange.f90->sourcefile~procedures_boxes_volume_exchange_factory.f90 sourcefile~classes_volume_change_method.f90->sourcefile~types_markov_chain_explorer_wrapper.f90 sourcefile~classes_volume_change_method.f90->sourcefile~procedures_exploration_inquirers.f90 sourcefile~procedures_volume_change_method_factory.f90 procedures_volume_change_method_factory.f90 sourcefile~classes_volume_change_method.f90->sourcefile~procedures_volume_change_method_factory.f90 sourcefile~procedures_box_volume_change_factory.f90 procedures_box_volume_change_factory.f90 sourcefile~classes_box_volume_change.f90->sourcefile~procedures_box_volume_change_factory.f90 sourcefile~procedures_generating_algorithms_factory.f90 procedures_generating_algorithms_factory.f90 sourcefile~procedures_boxes_volume_exchange_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_generating_algorithms_factory.f90->sourcefile~plmc_generate.f90 sourcefile~procedures_markov_chain_generator_factory.f90 procedures_markov_chain_generator_factory.f90 sourcefile~procedures_generating_algorithms_factory.f90->sourcefile~procedures_markov_chain_generator_factory.f90 sourcefile~procedures_volume_change_method_factory.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_box_volume_change_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_maximum_boxes_compression_explorer_factory.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_plmc_resetter.f90->sourcefile~plmc_explore.f90 sourcefile~procedures_plmc_resetter.f90->sourcefile~plmc_generate.f90 sourcefile~procedures_cells_memento.f90->sourcefile~classes_boxes_volume_exchange.f90 sourcefile~procedures_cells_memento.f90->sourcefile~classes_volume_change_method.f90 sourcefile~procedures_cells_memento.f90->sourcefile~classes_box_volume_change.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~procedures_short_interactions_factory.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_boxes_volume_exchange.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_volume_change_method.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_box_volume_change.f90 sourcefile~classes_box_particle_exchange.f90 classes_box_particle_exchange.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_box_particle_exchange.f90 sourcefile~classes_boxes_particles_swap.f90 classes_boxes_particles_swap.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_boxes_particles_swap.f90 sourcefile~classes_box_particles_swap.f90 classes_box_particles_swap.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_box_particles_swap.f90 sourcefile~classes_particle_insertion_method.f90 classes_particle_insertion_method.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_particle_insertion_method.f90 sourcefile~types_physical_model_wrapper.f90 types_physical_model_wrapper.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~types_physical_model_wrapper.f90 sourcefile~classes_box_particle_move.f90 classes_box_particle_move.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_box_particle_move.f90 sourcefile~classes_boxes_particle_teleportation.f90 classes_boxes_particle_teleportation.f90 sourcefile~types_short_interactions_wrapper.f90->sourcefile~classes_boxes_particle_teleportation.f90 sourcefile~procedures_box_particle_exchange_factory.f90 procedures_box_particle_exchange_factory.f90 sourcefile~classes_box_particle_exchange.f90->sourcefile~procedures_box_particle_exchange_factory.f90 sourcefile~procedures_boxes_particles_swap_factory.f90 procedures_boxes_particles_swap_factory.f90 sourcefile~classes_boxes_particles_swap.f90->sourcefile~procedures_boxes_particles_swap_factory.f90 sourcefile~procedures_box_particles_swap_factory.f90 procedures_box_particles_swap_factory.f90 sourcefile~classes_box_particles_swap.f90->sourcefile~procedures_box_particles_swap_factory.f90 sourcefile~classes_particle_insertion_method.f90->sourcefile~types_markov_chain_explorer_wrapper.f90 sourcefile~classes_particle_insertion_method.f90->sourcefile~procedures_exploration_inquirers.f90 sourcefile~procedures_particle_insertion_method_factory.f90 procedures_particle_insertion_method_factory.f90 sourcefile~classes_particle_insertion_method.f90->sourcefile~procedures_particle_insertion_method_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_dipolar_neighbourhoods_visitors_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~plmc_explore.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_physical_model_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_plmc_visitor.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~plmc_generate.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_boxes_volume_exchange_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_markov_chain_generator_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_volume_change_method_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_box_volume_change_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_maximum_boxes_compression_explorer_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_plmc_resetter.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_box_particle_exchange_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_boxes_particles_swap_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_box_particles_swap_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_particle_insertion_method_factory.f90 sourcefile~procedures_boxes_particle_teleportation_factory.f90 procedures_boxes_particle_teleportation_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_boxes_particle_teleportation_factory.f90 sourcefile~procedures_box_particle_move_factory.f90 procedures_box_particle_move_factory.f90 sourcefile~types_physical_model_wrapper.f90->sourcefile~procedures_box_particle_move_factory.f90 sourcefile~classes_box_particle_move.f90->sourcefile~procedures_box_particle_move_factory.f90 sourcefile~classes_boxes_particle_teleportation.f90->sourcefile~procedures_boxes_particle_teleportation_factory.f90 sourcefile~procedures_box_particle_exchange_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_boxes_particles_swap_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_box_particles_swap_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_particle_insertion_method_factory.f90->sourcefile~procedures_markov_chain_explorer_factory.f90 sourcefile~procedures_boxes_particle_teleportation_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_box_particle_move_factory.f90->sourcefile~procedures_generating_algorithms_factory.f90 sourcefile~procedures_exchange_visitors.f90->sourcefile~classes_box_particle_exchange.f90 sourcefile~procedures_exchange_visitors.f90->sourcefile~classes_boxes_particle_teleportation.f90 sourcefile~procedures_transmutation_updaters.f90->sourcefile~classes_boxes_particles_swap.f90 sourcefile~procedures_transmutation_updaters.f90->sourcefile~classes_box_particles_swap.f90 sourcefile~procedures_exchange_updaters.f90->sourcefile~classes_box_particle_exchange.f90 sourcefile~procedures_exchange_updaters.f90->sourcefile~classes_boxes_particle_teleportation.f90 sourcefile~procedures_transmutation_visitors.f90->sourcefile~classes_boxes_particles_swap.f90 sourcefile~procedures_transmutation_visitors.f90->sourcefile~classes_box_particles_swap.f90
Help


Source Code

module classes_neighbour_cells

use, intrinsic :: iso_fortran_env, only: DP => REAL64
use data_constants, only: num_dimensions, real_zero
use data_cells, only: nums_local_cells
use procedures_errors, only: error_exit
use procedures_checks, only: check_positive
use classes_parallelepiped_domain, only: Abstract_Parallelepiped_Domain
use classes_pair_potential, only: Abstract_Pair_Potential
use classes_hard_contact, only: Abstract_Hard_Contact
use classes_dipoles_neighbourhood, only: Abstract_Dipolar_Neighbourhood

implicit none

private

    type, abstract, public :: Abstract_Neighbour_Cells
    private
        class(Abstract_Parallelepiped_Domain), pointer :: accessible_domain => null()
        real(DP) :: max_distance = 0._DP
        integer :: nums(num_dimensions) = 0
        real(DP) :: size(num_dimensions) = 0._DP
        integer, dimension(num_dimensions) :: global_lbounds = 0, global_ubounds = 0
        logical :: skip_bottom_layer(-nums_local_cells(3)/2:nums_local_cells(3)/2) = .false.
        logical :: skip_top_layer(-nums_local_cells(3)/2:nums_local_cells(3)/2) = .false.
        integer, allocatable :: neighbours(:, :, :, :, :, :, :)
    contains
        procedure :: construct => Abstract_construct
        procedure :: destroy => Abstract_destroy
        procedure :: resize_only => Abstract_resize_only
        procedure :: reset => Abstract_reset
        procedure :: get_global_lbounds => Abstract_get_global_lbounds
        procedure :: get_global_ubounds => Abstract_get_global_ubounds
        procedure :: get => Abstract_get
        procedure :: skip => Abstract_skip
        procedure :: is_inside => Abstract_is_inside
        procedure :: index => Abstract_index
        procedure(Abstract_set_skip_layers), private, deferred :: set_skip_layers
        procedure(Abstract_check_nums), private, deferred :: check_nums
        procedure, private :: check_size => Abstract_check_size
        procedure, private :: set_neighbours => Abstract_set_neighbours
    end type Abstract_Neighbour_Cells

    abstract interface

        pure subroutine Abstract_set_skip_layers(this)
        import :: Abstract_Neighbour_Cells
            class(Abstract_Neighbour_Cells), intent(inout) :: this
        end subroutine Abstract_set_skip_layers

        subroutine Abstract_check_nums(this)
        import :: Abstract_Neighbour_Cells
            class(Abstract_Neighbour_Cells), intent(in) :: this
        end subroutine Abstract_check_nums

    end interface

    type, extends(Abstract_Neighbour_Cells), public :: XYZ_PBC_Neighbour_Cells
    contains
        procedure, private :: set_skip_layers => XYZ_set_skip_layers
        procedure, private :: check_nums => XYZ_check_nums
    end type XYZ_PBC_Neighbour_Cells

    type, extends(Abstract_Neighbour_Cells), public :: XY_PBC_Neighbour_Cells
    contains
        procedure, private :: set_skip_layers => XY_set_skip_layers
        procedure, private :: check_nums => XY_check_nums
    end type XY_PBC_Neighbour_Cells

    type, extends(Abstract_Neighbour_Cells), public :: Null_Neighbour_Cells
    contains
        procedure :: construct => Null_construct
        procedure :: destroy => Null_destroy
        procedure :: resize_only => Null_resize_only
        procedure :: reset => Null_reset
        procedure :: get_global_lbounds => Null_get_global_bounds
        procedure :: get_global_ubounds => Null_get_global_bounds
        procedure :: get => Null_get
        procedure :: skip => Null_skip
        procedure :: is_inside => Null_is_inside
        procedure :: index => Null_index
        procedure, private :: check_nums => Null_check_nums
        procedure, private :: set_skip_layers => Null_set_skip_layers
    end type Null_Neighbour_Cells

    type, public :: Neighbour_Cells_Wrapper
        class(Abstract_Neighbour_Cells), allocatable :: cells
    end type Neighbour_Cells_Wrapper

    type, public :: Neighbour_Cells_Line
        type(Neighbour_Cells_Wrapper), allocatable :: line(:)
    end type Neighbour_Cells_Line

contains

!implementation Abstract_Neighbour_Cells

    subroutine Abstract_construct(this, accessible_domain, pair_potential, hard_contact, &
        dipolar_neighbourhood)
        class(Abstract_Neighbour_Cells), intent(out) :: this
        class(Abstract_Parallelepiped_Domain), target, intent(in) :: accessible_domain
        class(Abstract_Pair_Potential), intent(in) :: pair_potential
        class(Abstract_Hard_Contact), intent(in) :: hard_contact
        class(Abstract_Dipolar_Neighbourhood), intent(in) :: dipolar_neighbourhood

        this%accessible_domain => accessible_domain
        this%max_distance = pair_potential%get_max_distance()
        if (this%max_distance - pair_potential%get_min_distance() < hard_contact%get_width()) then
            this%max_distance = this%max_distance + hard_contact%get_width()
        end if
        if (this%max_distance < dipolar_neighbourhood%get_max_distance()) then
            this%max_distance = dipolar_neighbourhood%get_max_distance()
        end if
        call check_positive("Abstract_Neighbour_Cells: construct", "this%max_distance", this%&
            max_distance)
        call this%set_skip_layers()
    end subroutine Abstract_construct

    subroutine Abstract_destroy(this)
        class(Abstract_Neighbour_Cells), intent(inout) :: this

        if (allocated(this%neighbours)) deallocate(this%neighbours)
        this%accessible_domain => null()
    end subroutine Abstract_destroy

    pure logical function Abstract_resize_only(this) result(resize_only)
        class(Abstract_Neighbour_Cells), intent(in) :: this

        resize_only = all(this%nums == floor(this%accessible_domain%get_size() / this%max_distance))
    end function Abstract_resize_only

    subroutine Abstract_reset(this)
        class(Abstract_Neighbour_Cells), intent(inout) :: this

        integer :: nums(num_dimensions)

        nums = floor(this%accessible_domain%get_size() / this%max_distance)
        if (all(this%nums == nums)) then
            this%size = this%accessible_domain%get_size() / real(this%nums, DP)
            call this%check_size()
            return
        end if
        this%nums = nums
        call this%check_nums()
        this%size = this%accessible_domain%get_size() / real(this%nums, DP)
        call this%check_size()

        this%global_lbounds = -this%nums/2
        this%global_ubounds = this%global_lbounds + this%nums - 1

        if (allocated(this%neighbours)) deallocate(this%neighbours)
        allocate(this%neighbours(3, -nums_local_cells(1)/2:nums_local_cells(1)/2, &
                                    -nums_local_cells(2)/2:nums_local_cells(2)/2, &
                                    -nums_local_cells(3)/2:nums_local_cells(3)/2, &
                                    this%global_lbounds(1):this%global_ubounds(1), &
                                    this%global_lbounds(2):this%global_ubounds(2), &
                                    this%global_lbounds(3):this%global_ubounds(3)))
        call this%set_neighbours()
    end subroutine Abstract_reset

    subroutine Abstract_check_size(this)
        class(Abstract_Neighbour_Cells), intent(in) :: this

        real(DP) :: box_modulo_cell(num_dimensions)

        box_modulo_cell = modulo(this%accessible_domain%get_size(), this%size)
        if (any(box_modulo_cell > real_zero .and. abs(box_modulo_cell-this%size) > real_zero)) then
            call error_exit("Abstract_Neighbour_Cells: check_size: "//&
                            "this%size is not a divisor of accessible_domain%get_size()")
        end if
        ! rewrite to be smarter: handle exception
    end subroutine Abstract_check_size

    subroutine Abstract_set_neighbours(this)
        class(Abstract_Neighbour_Cells), intent(inout) :: this

        integer :: global_i1, global_i2, global_i3
        logical :: at_bottom_layer, at_top_layer
        integer :: local_i1, local_i2, local_i3
        integer :: ijk_cell(num_dimensions)

        this%neighbours = 0
        do global_i3 = this%global_lbounds(3), this%global_ubounds(3)
            at_bottom_layer = (global_i3 == this%global_lbounds(3))
            at_top_layer = (global_i3 == this%global_ubounds(3))
        do global_i2 = this%global_lbounds(2), this%global_ubounds(2)
        do global_i1 = this%global_lbounds(1), this%global_ubounds(1)
            do local_i3 = -nums_local_cells(3)/2, nums_local_cells(3)/2
                if (this%skip(at_bottom_layer, at_top_layer, local_i3)) cycle
            do local_i2 = -nums_local_cells(2)/2, nums_local_cells(2)/2
            do local_i1 = -nums_local_cells(1)/2, nums_local_cells(1)/2
                ijk_cell = [global_i1, global_i2, global_i3] + [local_i1, local_i2, local_i3]
                ijk_cell = pbc_3d_index(ijk_cell, this%nums)
                this%neighbours(:, local_i1, local_i2, local_i3, global_i1, global_i2, global_i3) =&
                    ijk_cell
            end do
            end do
            end do
        end do
        end do
        end do
    end subroutine Abstract_set_neighbours

    pure function Abstract_get_global_lbounds(this) result(lbounds)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        integer :: lbounds(num_dimensions)

        lbounds = this%global_lbounds
    end function Abstract_get_global_lbounds

    pure function Abstract_get_global_ubounds(this) result(ubounds)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        integer :: ubounds(num_dimensions)

        ubounds = this%global_ubounds
    end function Abstract_get_global_ubounds

    pure function Abstract_get(this, local_i1, local_i2, local_i3, global_i1, global_i2, global_i3)&
        result(neighbour)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        integer, intent(in) :: local_i1, local_i2, local_i3
        integer, intent(in) :: global_i1, global_i2, global_i3
        integer :: neighbour(num_dimensions)

        neighbour = this%neighbours(:, local_i1, local_i2, local_i3, global_i1, global_i2, &
            global_i3)
    end function Abstract_get

    pure logical function Abstract_skip(this, at_bottom_layer, at_top_layer, local_i3) result(skip)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        logical, intent(in) :: at_bottom_layer, at_top_layer
        integer, intent(in) :: local_i3

        skip = (at_bottom_layer .and. this%skip_bottom_layer(local_i3)) .or. &
            (at_top_layer .and. this%skip_top_layer(local_i3))
    end function Abstract_skip

    pure logical function Abstract_is_inside(this, position) result(is_inside)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        real(DP), intent(in) :: position(:)

        is_inside = this%accessible_domain%is_inside(position)
    end function Abstract_is_inside

    pure function Abstract_index(this, position) result(index)
        class(Abstract_Neighbour_Cells), intent(in) :: this
        real(DP), intent(in) :: position(:)
        integer :: index(num_dimensions)

        where (mod(this%nums, 2) == 0)
            index = floor(position/this%size)
        elsewhere
            index = nint(position/this%size)
        end where
    end function Abstract_index

!end implementation Abstract_Neighbour_Cells

    !> @note Periodic index, found heuristically.
    pure function pbc_3d_index(ijk_cell, nums_cells)
        integer, intent(in) :: ijk_cell(:), nums_cells(:)
        integer :: pbc_3d_index(3)

        pbc_3d_index = modulo(ijk_cell + nums_cells/2, nums_cells) - nums_cells/2
    end function pbc_3d_index

!implementation XYZ_PBC_Neighbour_Cells

    pure subroutine XYZ_set_skip_layers(this)
        class(XYZ_PBC_Neighbour_Cells), intent(inout) :: this

        this%skip_top_layer = .false.
        this%skip_bottom_layer = .false.
    end subroutine XYZ_set_skip_layers

    subroutine XYZ_check_nums(this)
        class(XYZ_PBC_Neighbour_Cells), intent(in) :: this

        if (any(this%nums < nums_local_cells)) then
            call error_exit("XYZ_PBC_Neighbour_Cells: this%nums is too small.")
        end if
    end subroutine XYZ_check_nums

!end implementation XYZ_PBC_Neighbour_Cells

!implementation XY_PBC_Neighbour_Cells

    pure subroutine XY_set_skip_layers(this)
        class(XY_PBC_Neighbour_Cells), intent(inout) :: this

        integer :: local_i3

        do local_i3 = -nums_local_cells(3)/2, nums_local_cells(3)/2
            if (local_i3 == 1) then
                this%skip_top_layer(local_i3) = .true.
            else
                this%skip_top_layer(local_i3) = .false.
            end if
            if (local_i3 == -1) then
                this%skip_bottom_layer(local_i3) = .true.
            else
                this%skip_bottom_layer(local_i3) = .false.
            end if
        end do
    end subroutine XY_set_skip_layers

    subroutine XY_check_nums(this)
        class(XY_PBC_Neighbour_Cells), intent(in) :: this

        if (any(this%nums(1:2) < nums_local_cells(1:2)) .or. this%nums(3) < 1) then
            call error_exit("XY_PBC_Neighbour_Cells: this%nums is too small.")
        end if
    end subroutine XY_check_nums

!end implementation XY_PBC_Neighbour_Cells

!implementation Null_Neighbour_Cells

    subroutine Null_construct(this, accessible_domain, pair_potential, hard_contact, &
        dipolar_neighbourhood)
        class(Null_Neighbour_Cells), intent(out) :: this
        class(Abstract_Parallelepiped_Domain), target, intent(in) :: accessible_domain
        class(Abstract_Pair_Potential), intent(in) :: pair_potential
        class(Abstract_Hard_Contact), intent(in) :: hard_contact
        class(Abstract_Dipolar_Neighbourhood), intent(in) :: dipolar_neighbourhood
    end subroutine Null_construct

    subroutine Null_destroy(this)
        class(Null_Neighbour_Cells), intent(inout) :: this
    end subroutine Null_destroy

    pure logical function Null_resize_only(this) result(resize_only)
        class(Null_Neighbour_Cells), intent(in) :: this
        resize_only = .false.
    end function Null_resize_only

    subroutine Null_reset(this)
        class(Null_Neighbour_Cells), intent(inout) :: this
    end subroutine Null_reset

    pure function Null_get_global_bounds(this) result(bounds)
        class(Null_Neighbour_Cells), intent(in) :: this
        integer :: bounds(num_dimensions)
        bounds = 0
    end function Null_get_global_bounds

    pure function Null_get(this, local_i1, local_i2, local_i3, global_i1, global_i2, global_i3)&
        result(neighbour)
        class(Null_Neighbour_Cells), intent(in) :: this
        integer, intent(in) :: local_i1, local_i2, local_i3
        integer, intent(in) :: global_i1, global_i2, global_i3
        integer :: neighbour(num_dimensions)
        neighbour = 0
    end function Null_get

    pure logical function Null_skip(this, at_bottom_layer, at_top_layer, local_i3) result(skip)
        class(Null_Neighbour_Cells), intent(in) :: this
        logical, intent(in) :: at_bottom_layer, at_top_layer
        integer, intent(in) :: local_i3
        skip = .false.
    end function Null_skip

    pure logical function Null_is_inside(this, position) result(is_inside)
        class(Null_Neighbour_Cells), intent(in) :: this
        real(DP), intent(in) :: position(:)
        is_inside = .false.
    end function Null_is_inside

    pure function Null_index(this, position) result(index)
        class(Null_Neighbour_Cells), intent(in) :: this
        real(DP), intent(in) :: position(:)
        integer :: index(num_dimensions)
        index = 0
    end function Null_index

    pure subroutine Null_set_skip_layers(this)
        class(Null_Neighbour_Cells), intent(inout) :: this
    end subroutine Null_set_skip_layers

    subroutine Null_check_nums(this)
        class(Null_Neighbour_Cells), intent(in) :: this
    end subroutine Null_check_nums

!end implementation Null_Neighbour_Cells

end module classes_neighbour_cells

classes_average_num_particles.f90 classes_beta_pressure.f90 classes_beta_pressure_excess.f90 classes_box_particle_exchange.f90 classes_box_particle_move.f90 classes_box_particles_swap.f90 classes_box_size_checker.f90 classes_box_size_memento.f90 classes_box_volume_change.f90 classes_boxes_particle_teleportation.f90 classes_boxes_particles_swap.f90 classes_boxes_volume_exchange.f90 classes_changed_box_size.f90 classes_changed_box_size_ratio.f90 classes_changes_success_writer.f90 classes_complete_coordinates_reader.f90 classes_complete_coordinates_writer.f90 classes_component_chemical_potential.f90 classes_component_coordinates.f90 classes_component_coordinates_reader.f90 classes_component_coordinates_writer.f90 classes_component_dipole_moments.f90 classes_coordinates.f90 classes_coordinates_copier.f90 classes_density_explorer.f90 classes_des_convergence_parameter.f90 classes_des_real_component.f90 classes_des_real_pair.f90 classes_des_reci_structure.f90 classes_des_reci_visitor.f90 classes_des_reci_weight.f90 classes_des_self_component.f90 classes_des_surf_mixture.f90 classes_dipolar_interactions_facade.f90 classes_dipolar_neighbourhood.f90 classes_dipolar_neighbourhoods_visitor.f90 classes_dirac_distribution_plus.f90 classes_directed_graph_writer.f90 classes_dlc_structures.f90 classes_dlc_visitor.f90 classes_dlc_weight.f90 classes_exchanged_boxes_size.f90 classes_external_field.f90 classes_field_expression.f90 classes_floor_penetration.f90 classes_generating_algorithm.f90 classes_hard_contact.f90 classes_hetero_couples.f90 classes_line_writer.f90 classes_maximum_box_compression.f90 classes_maximum_box_compression_explorer.f90 classes_min_distance.f90 classes_mixture_total_moment.f90 classes_move_tuner.f90 classes_moved_coordinates.f90 classes_neighbour_cells.f90 classes_num_particles.f90 classes_number_to_string.f90 classes_pair_potential.f90 classes_parallelepiped_domain.f90 classes_particle_insertion_method.f90 classes_periodic_box.f90 classes_permittivity.f90 classes_plmc_propagator.f90 classes_potential_expression.f90 classes_radial_explorer.f90 classes_random_coordinates.f90 classes_random_orientation.f90 classes_random_position.f90 classes_real_writer.f90 classes_reciprocal_lattice.f90 classes_rectangle_writer.f90 classes_rotated_orientations.f90 classes_short_pairs_visitor.f90 classes_structure_factor.f90 classes_structure_visitor.f90 classes_temperature.f90 classes_tower_sampler.f90 classes_translated_positions.f90 classes_triangle_writer.f90 classes_tunable_move.f90 classes_visitable_cells.f90 classes_visitable_cells_memento.f90 classes_visitable_list.f90 classes_visitable_walls.f90 classes_volume_change_method.f90 classes_walls_visitor.f90 data_cells.f90 data_constants.f90 data_input_prefixes.f90 data_output_objects.f90 data_strings.f90 density.f90 module_changes_success.f90 module_list_node.f90 module_move_tuning.f90 plmc_explore.f90 plmc_generate.f90 procedures_average_nums_particles_factory.f90 procedures_beta_pressure_factory.f90 procedures_beta_pressures_excess_factory.f90 procedures_box_particle_exchange_factory.f90 procedures_box_particle_move_factory.f90 procedures_box_particles_swap_factory.f90 procedures_box_size.f90 procedures_box_size_memento_factory.f90 procedures_box_volume_change_factory.f90 procedures_boxes_particle_teleportation_factory.f90 procedures_boxes_particles_swap_factory.f90 procedures_boxes_size_checker_factory.f90 procedures_boxes_volume_exchange_factory.f90 procedures_cells_memento.f90 procedures_centered_block_micro.f90 procedures_changed_boxes_size_factory.f90 procedures_changed_boxes_size_ratio_factory.f90 procedures_changes_component_factory.f90 procedures_changes_factory.f90 procedures_changes_properties.f90 procedures_changes_success_writer_factory.f90 procedures_checks.f90 procedures_command_arguments.f90 procedures_complete_coordinates_reader.f90 procedures_complete_coordinates_reader_factory.f90 procedures_complete_coordinates_writer_factory.f90 procedures_component_chemical_potential_factory.f90 procedures_component_coordinates_factory.f90 procedures_component_coordinates_reader_factory.f90 procedures_component_coordinates_writer_factory.f90 procedures_component_dipole_moments_factory.f90 procedures_component_factory.f90 procedures_coordinates_copier_factory.f90 procedures_coordinates_micro.f90 procedures_coordinates_reader.f90 procedures_density_explorer_factory.f90 procedures_des_convergence_parameter_factory.f90 procedures_des_real_component_factory.f90 procedures_des_real_pair_factory.f90 procedures_des_reci_structure_factory.f90 procedures_des_reci_visitor_factory.f90 procedures_des_reci_weight_factory.f90 procedures_des_surf_mixture_factory.f90 procedures_dipolar_interactions_facades_factory.f90 procedures_dipolar_interactions_factory.f90 procedures_dipolar_interactions_micro.f90 procedures_dipolar_interactions_resetter.f90 procedures_dipolar_interactions_visitor.f90 procedures_dipolar_neighbourhoods_factory.f90 procedures_dipolar_neighbourhoods_visitors_factory.f90 procedures_dipoles_field_interaction.f90 procedures_dirac_distribution_plus.f90 procedures_directed_graph_writer_factory.f90 procedures_dlc_structures_factory.f90 procedures_dlc_visitor_factory.f90 procedures_dlc_weight_factory.f90 procedures_elementary_geometry.f90 procedures_elementary_statistics.f90 procedures_energies_writers_factory.f90 procedures_environment_factory.f90 procedures_environment_inquirers.f90 procedures_errors.f90 procedures_exchange_updaters.f90 procedures_exchange_visitors.f90 procedures_exchanged_boxes_size_factory.f90 procedures_exploration_inquirers.f90 procedures_exploring_observables_factory.f90 procedures_exploring_writers_factory.f90 procedures_external_fields_factory.f90 procedures_field_expression_factory.f90 procedures_field_expression_micro.f90 procedures_floor_penetration_factory.f90 procedures_generating_algorithms_factory.f90 procedures_generating_observables_factory.f90 procedures_generating_writers_factory.f90 procedures_hard_contact_factory.f90 procedures_hetero_couples_factory.f90 procedures_json_data_factory.f90 procedures_json_reports_factory.f90 procedures_line_writer_factory.f90 procedures_logical_factory.f90 procedures_markov_chain_explorer_factory.f90 procedures_markov_chain_generator_factory.f90 procedures_maximum_box_compression_factory.f90 procedures_maximum_boxes_compression_explorer_factory.f90 procedures_metropolis_algorithm.f90 procedures_min_distance_factory.f90 procedures_mixture_inquirers.f90 procedures_mixture_properties.f90 procedures_mixture_total_moments_factory.f90 procedures_move_tuner_factory.f90 procedures_moved_coordinates_factory.f90 procedures_neighbour_cells_factory.f90 procedures_num_particles_factory.f90 procedures_observables_changes_factory.f90 procedures_observables_energies_factory.f90 procedures_observables_factory.f90 procedures_pair_potential_factory.f90 procedures_parallelepiped_domain_macro.f90 procedures_parallelepiped_domain_micro.f90 procedures_parallelepiped_domains_factory.f90 procedures_particle_insertion_method_factory.f90 procedures_periodic_boxes_factory.f90 procedures_permittivity_factory.f90 procedures_physical_model_factory.f90 procedures_plmc_help.f90 procedures_plmc_iterations.f90 procedures_plmc_propagator_factory.f90 procedures_plmc_resetter.f90 procedures_plmc_visitor.f90 procedures_plmc_writer.f90 procedures_potential_expression_factory.f90 procedures_property_inquirers.f90 procedures_radial_explorer_factory.f90 procedures_random_coordinates_factory.f90 procedures_random_number.f90 procedures_random_seed_factory.f90 procedures_readers_factory.f90 procedures_real_writer_factory.f90 procedures_reals_factory.f90 procedures_reciprocal_lattices_factory.f90 procedures_rectangle_writer_factory.f90 procedures_selectors_resetters.f90 procedures_short_interactions_factory.f90 procedures_short_interactions_inquirers.f90 procedures_short_interactions_resetter.f90 procedures_short_interactions_visitor.f90 procedures_short_pairs_visitors_factory.f90 procedures_string_factory.f90 procedures_temperature_factory.f90 procedures_tower_sampler_factory.f90 procedures_transmutation_updaters.f90 procedures_transmutation_visitors.f90 procedures_triangle_observables.f90 procedures_triangle_writer_factory.f90 procedures_visit_condition.f90 procedures_visitable_cells_factory.f90 procedures_visitable_cells_memento_factory.f90 procedures_visitable_list_factory.f90 procedures_visitable_walls_factory.f90 procedures_volume_change_method_factory.f90 procedures_walls_visitors_factory.f90 procedures_writers_inquirers.f90 radial.f90 types_cells_wrapper.f90 types_changes_component_wrapper.f90 types_changes_success_writer_selector.f90 types_changes_wrapper.f90 types_component_coordinates_reader_selector.f90 types_component_coordinates_writer_selector.f90 types_component_wrapper.f90 types_dipolar_interactions_dynamic_wrapper.f90 types_dipolar_interactions_static_wrapper.f90 types_energies_writers.f90 types_environment_wrapper.f90 types_exploring_io.f90 types_exploring_writers_wrapper.f90 types_generating_io.f90 types_generating_observables_wrapper.f90 types_generating_writers_wrapper.f90 types_json_report.f90 types_logical_wrapper.f90 types_markov_chain_explorer_wrapper.f90 types_markov_chain_generator_wrapper.f90 types_mixture_wrapper.f90 types_move_tuner_parameters.f90 types_observables_changes.f90 types_observables_energies.f90 types_particle_wrapper.f90 types_physical_model_wrapper.f90 types_potential_domain.f90 types_potential_domain_selector.f90 types_raw_coordinates.f90 types_readers_wrapper.f90 types_real_wrapper.f90 types_short_interactions_wrapper.f90 types_string_wrapper.f90