![]() |
ESPResSo 3.2.0-11-g9950804-git
Extensible Simulation Package for Soft Matter Research
|
Domain decomposition for parallel computing. More...
#include <mpi.h>#include <stdio.h>#include <stdlib.h>#include <string.h>#include <math.h>#include "utils.h"#include "grid.h"#include "communication.h"#include "verlet.h"#include "cells.h"#include "interaction_data.h"
Go to the source code of this file.
Defines | |
| #define | MAX_INTERACTION_RANGE 1e100 |
Functions | |
| void | init_node_grid () |
| Make sure that the node grid is set, eventually determine one automatically. | |
| int | node_grid_is_set () |
| return wether node grid was set. | |
| int | map_position_node_array (double pos[3]) |
| map a spatial position to the node grid | |
| void | calc_node_neighbors (int node) |
| fill neighbor lists of node. | |
| void | grid_changed_box_l () |
| called from mpi_bcast_parameter . | |
| void | grid_changed_n_nodes () |
| called from mpi_bcast_parameter . | |
| void | calc_minimal_box_dimensions () |
| Calculates the smallest box and local box dimensions for periodic directions. | |
| void | calc_2d_grid (int n, int grid[3]) |
| calculate most square 2d grid. | |
| int | map_3don2d_grid (int g3d[3], int g2d[3], int mult[3]) |
| calculate 'best' mapping between a 2d and 3d grid. | |
| void | rescale_boxl (int dir, double d_new) |
| rescales the box in dimension 'dir' to the new value 'd_new', and rescales the particles accordingly | |
Variables | |
| int | node_grid [3] = { 0, 0, 0} |
| The number of nodes in each spatial dimension. | |
| int | node_pos [3] = {-1,-1,-1} |
| position of node in node grid | |
| int | node_neighbors [6] = {0, 0, 0, 0, 0, 0} |
| the six nearest neighbors of a node in the node grid. | |
| int | boundary [6] = {0, 0, 0, 0, 0, 0} |
| where to fold particles that leave local box in direction i. | |
| int | periodic = 7 |
| Flags for all three dimensions wether pbc are applied (default). | |
| double | box_l [3] = {1, 1, 1} |
| Simulation box dimensions. | |
| double | box_l_i [3] = {1, 1, 1} |
| 1 / box dimensions. | |
| double | min_box_l |
| Smallest simulation box dimension (box_l). | |
| double | local_box_l [3] = {1, 1, 1} |
| Dimensions of the box a single node is responsible for. | |
| double | min_local_box_l |
| Smallest local simulation box dimension (local_box_l). | |
| double | my_left [3] = {0, 0, 0} |
| Left (bottom, front) corner of this nodes local box. | |
| double | my_right [3] = {1, 1, 1} |
| Right (top, back) corner of this nodes local box. | |
Domain decomposition for parallel computing.
For more information on the domain decomposition, see grid.h.
Definition in file grid.c.
| #define MAX_INTERACTION_RANGE 1e100 |
Definition at line 42 of file grid.c.
Referenced by calc_minimal_box_dimensions().
| void calc_2d_grid | ( | int | n, |
| int | grid[3] | ||
| ) |
calculate most square 2d grid.
Definition at line 186 of file grid.c.
References i.
Referenced by dfft_init(), and fft_init().
| void calc_minimal_box_dimensions | ( | ) |
Calculates the smallest box and local box dimensions for periodic directions.
This is needed to check if the interaction ranges are compatible with the box dimensions and the node grid. see also box_l, local_box_l, min_box_l and min_local_box_l. Remark: In the apreiodic case min_box_l is set to 2 * MAX_INTERACTION_RANGE .
Definition at line 175 of file grid.c.
References box_l, dmin(), i, local_box_l, MAX_INTERACTION_RANGE, min_box_l, and min_local_box_l.
Referenced by grid_changed_box_l().
| void calc_node_neighbors | ( | int | node | ) |
fill neighbor lists of node.
Calculates the numbers of the 6 nearest neighbors for a node and stors them in node_neighbors.
| node | number of the node. |
Definition at line 96 of file grid.c.
References boundary, comm_cart, GRID_TRACE, map_node_array(), MPI_Cart_shift(), node_grid, node_neighbors, node_pos, and this_node.
Referenced by grid_changed_n_nodes().
| void grid_changed_box_l | ( | ) |
called from mpi_bcast_parameter .
Definition at line 124 of file grid.c.
References box_l, box_l_i, calc_minimal_box_dimensions(), GRID_TRACE, i, local_box_l, my_left, my_right, node_grid, node_pos, and this_node.
Referenced by grid_changed_n_nodes(), on_parameter_change(), and propagate_press_box_pos_and_rescale_npt().
| void grid_changed_n_nodes | ( | ) |
called from mpi_bcast_parameter .
Definition at line 148 of file grid.c.
References boundary, calc_node_neighbors(), comm_cart, grid_changed_box_l(), GRID_TRACE, MPI_Cart_coords(), MPI_Cart_create(), MPI_Comm_free(), MPI_Comm_rank(), MPI_COMM_WORLD, node_grid, node_neighbors, node_pos, and this_node.
Referenced by init_node_grid(), and on_parameter_change().
| void init_node_grid | ( | ) |
Make sure that the node grid is set, eventually determine one automatically.
Definition at line 64 of file grid.c.
References CELL_FLAG_GRIDCHANGED, cells_on_geometry_change(), and grid_changed_n_nodes().
Referenced by on_program_start().
| int map_3don2d_grid | ( | int | g3d[3], |
| int | g2d[3], | ||
| int | mult[3] | ||
| ) |
calculate 'best' mapping between a 2d and 3d grid.
This we need for the communication from 3d domain decomposition to 2d row decomposition. The dimensions of the 2d grid are resorted, if necessary, in a way that they are multiples of the 3d grid dimensions.
| g3d | 3d grid. |
| g2d | 2d grid. |
| mult | factors between 3d and 2d grid dimensions |
Definition at line 196 of file grid.c.
References i.
Referenced by dfft_init(), and fft_init().
| int map_position_node_array | ( | double | pos[3] | ) |
map a spatial position to the node grid
Definition at line 75 of file grid.c.
References box_l_i, fold_position(), i, map_array_node(), and node_grid.
Referenced by dd_topology_init(), layered_topology_init(), and nsq_topology_init().
| int node_grid_is_set | ( | ) |
| void rescale_boxl | ( | int | dir, |
| double | d_new | ||
| ) |
rescales the box in dimension 'dir' to the new value 'd_new', and rescales the particles accordingly
Definition at line 220 of file grid.c.
References box_l, FIELD_BOXL, mpi_bcast_parameter(), and mpi_rescale_particles().
Referenced by tclcommand_change_volume().
| int boundary[6] = {0, 0, 0, 0, 0, 0} |
where to fold particles that leave local box in direction i.
Definition at line 51 of file grid.c.
Referenced by calc_node_neighbors(), dd_append_particles(), dd_exchange_and_sort_particles(), dd_prepare_comm(), dd_save_position_to_cell(), dd_update_communicators_w_boxl(), grid_changed_n_nodes(), lb_bounce_back(), lb_calc_local_fields(), lb_collide_stream(), lb_lbfluid_get_interpolated_velocity(), lb_lbfluid_print_boundary(), lb_lbfluid_print_vtk_boundary(), and prepare_halo_communication().
| double box_l[3] = {1, 1, 1} |
Simulation box dimensions.
Definition at line 54 of file grid.c.
Referenced by add_dipole_force(), add_ljangle_pair_force(), add_mdlc_energy_corrections(), add_mdlc_force_corrections(), add_mmm1d_coulomb_pair_force(), add_mmm2d_coulomb_pair_force(), adress_dw_dir(), adress_wf_vector(), analyze_fold_molecules(), analyze_rdfchain(), bilayer_density_profile(), buf_mindist4(), build_verlet_lists(), build_verlet_lists_and_calc_verlet_ia(), build_verlet_lists_and_calc_verlet_ia_iccp3m(), calc_angledist_param(), calc_diffusion_profile(), calc_fluctuations(), calc_minimal_box_dimensions(), calc_mmm2d_copy_pair_energy(), calc_part_distribution(), calc_rdf(), calc_rdf_adress(), calc_rdf_av(), calc_rdf_intermol_av(), calc_scaling(), calc_scaling2(), calc_structurefactor(), calc_trap_force(), calc_wallrdfyz(), calculate_maze_dist(), check_particle_consistency(), check_particles(), counterionsC(), create_free_volume_grid(), dd_create_cell_grid(), dd_exchange_and_sort_particles(), dd_position_to_cell(), dd_prepare_comm(), dd_save_position_to_cell(), dd_update_communicators_w_boxl(), density_profile_av(), dipole_energy(), distribute_tensors(), distto(), does_line_go_through_cube(), dp3m_adaptive_tune(), dp3m_calc_influence_function_energy(), dp3m_calc_influence_function_force(), dp3m_calc_kspace_forces(), dp3m_compute_constants_energy_dipolar(), dp3m_get_accuracy(), dp3m_init_a_ai_cao_cut(), dp3m_mc_time(), dp3m_sanity_checks(), dp3m_sanity_checks_boxl(), dp3m_scaleby_box_l(), dp3m_set_params(), dp3m_set_tune_params(), ELC_set_params(), ELC_setup_constants(), ELC_tune(), fold_coordinate(), get_DLC_dipolar(), get_DLC_energy_dipolar(), get_lipid_orients(), get_mi_vector(), getlength(), grid_changed_box_l(), image_sum_b(), image_sum_t(), init_lattice(), init_local_particle_force(), integrate_ensemble_init(), layered_get_mi_vector(), layered_prepare_comm(), lb_calc_average_rho(), lb_calc_densprof(), lb_calc_velprof(), lb_lbfluid_get_interpolated_velocity_global(), lb_lbfluid_load_checkpoint(), lb_lbfluid_print_boundary(), lb_lbfluid_print_velocity(), lb_lbfluid_print_vtk_boundary(), lb_lbfluid_print_vtk_velocity(), lb_lbfluid_save_checkpoint(), lb_lbfluid_set_agrid(), ljangle_pair_energy(), local_stress_tensor_calc(), maggs_compute_dipole_correction(), maggs_init(), maggs_sanity_checks(), magnetic_dipolar_direct_sum_calculations(), mdlc_set_params(), mdlc_tune(), mindist(), mindist3(), mindist4(), mmm1d_coulomb_pair_energy(), MMM1D_init(), MMM1D_setup_constants(), mmm1d_tune(), MMM2D_add_far(), MMM2D_dielectric_layers_energy_contribution(), MMM2D_dielectric_layers_force_contribution(), MMM2D_setup_constants(), MMM2D_tune_near(), nemd_init(), observable_lb_radial_velocity_profile(), observable_lb_velocity_profile(), observable_structure_factor(), p3m_adaptive_tune(), p3m_calc_dipole_term(), p3m_calc_influence_function_force(), p3m_calc_kspace_forces(), p3m_calc_kspace_stress(), p3m_init_a_ai_cao_cut(), p3m_k_space_error(), p3m_mc_time(), p3m_perform_aliasing_sums_energy(), p3m_perform_aliasing_sums_force(), p3m_print(), p3m_real_space_error(), p3m_sanity_checks_boxl(), p3m_scaleby_box_l(), p3m_set_params(), p3m_set_tune_params(), polymerC(), pressure_calc(), propagate_press_box_pos_and_rescale_npt(), recalc_global_maximal_nonbonded_cutoff(), rescale_boxl(), saltC(), setup_P(), setup_PQ(), setup_Q(), tclcallback_box_l(), tclcommand_adress_parse_set(), tclcommand_analyze_fluid_parse_velprof(), tclcommand_analyze_parse_density_profile_av(), tclcommand_analyze_parse_diffusion_profile(), tclcommand_analyze_parse_structurefactor(), tclcommand_analyze_parse_Vkappa(), tclcommand_change_volume(), tclcommand_crosslink(), tclcommand_inter_coulomb_parse_mmm1d(), tclcommand_inter_parse_ljangle(), tclcommand_lbnode(), tclcommand_nemd_parse_and_print_viscosity(), tclcommand_parse_profile(), tclcommand_parse_radial_profile(), unfold_position(), update_mol_pos_particle(), and z_energy().
| double box_l_i[3] = {1, 1, 1} |
1 / box dimensions.
Definition at line 55 of file grid.c.
Referenced by add_ljangle_pair_force(), analyze_fold_molecules(), calc_surface_term(), dp3m_adaptive_tune(), dp3m_calc_kspace_forces(), dp3m_scaleby_box_l(), dp3m_set_params(), dp3m_set_tune_params(), fold_coordinate(), get_mi_vector(), grid_changed_box_l(), layered_get_mi_vector(), ljangle_pair_energy(), map_lattice_to_node(), map_position_node_array(), p3m_adaptive_tune(), p3m_calc_dipole_term(), p3m_scaleby_box_l(), p3m_set_params(), p3m_set_tune_params(), and rod_energy().
| double local_box_l[3] = {1, 1, 1} |
Dimensions of the box a single node is responsible for.
Definition at line 57 of file grid.c.
Referenced by calc_minimal_box_dimensions(), dd_create_cell_grid(), dd_on_geometry_change(), dd_position_to_cell(), dd_save_position_to_cell(), dp3m_mc_time(), dp3m_sanity_checks_boxl(), force_and_velocity_check(), grid_changed_box_l(), init_lattice(), layered_topology_init(), map_position_to_lattice(), p3m_mc_time(), and p3m_sanity_checks_boxl().
| double min_box_l |
Smallest simulation box dimension (box_l).
Remark: with PARTIAL_PERIODIC, only the periodic directions are taken into account!
Definition at line 56 of file grid.c.
Referenced by calc_minimal_box_dimensions(), dp3m_adaptive_tune(), dp3m_mc_time(), p3m_adaptive_tune(), p3m_mc_time(), tclcommand_analyze_parse_distribution(), and tclcommand_analyze_parse_rdf().
| double min_local_box_l |
Smallest local simulation box dimension (local_box_l).
Remark: with PARTIAL_PERIODIC, only the periodic directions are taken into account!
Definition at line 58 of file grid.c.
Referenced by calc_minimal_box_dimensions(), dp3m_adaptive_tune(), dp3m_mc_time(), p3m_adaptive_tune(), and p3m_mc_time().
| double my_left[3] = {0, 0, 0} |
Left (bottom, front) corner of this nodes local box.
Definition at line 59 of file grid.c.
Referenced by calc_particle_lattice_ia(), calc_viscous_force(), dd_append_particles(), dd_create_cell_grid(), dd_exchange_and_sort_particles(), dd_position_to_cell(), dd_save_position_to_cell(), dp3m_assign_dipole(), dp3m_calc_local_ca_mesh(), grid_changed_box_l(), layered_append_particles(), layered_exchange_and_sort_particles(), layered_position_to_cell(), maggs_setup_local_lattice(), map_position_to_lattice(), p3m_assign_charge(), p3m_calc_local_ca_mesh(), setup_P(), setup_PQ(), and setup_Q().
| double my_right[3] = {1, 1, 1} |
Right (top, back) corner of this nodes local box.
Definition at line 60 of file grid.c.
Referenced by calc_particle_lattice_ia(), dd_create_cell_grid(), dd_exchange_and_sort_particles(), dp3m_assign_dipole(), dp3m_calc_local_ca_mesh(), grid_changed_box_l(), layered_append_particles(), layered_exchange_and_sort_particles(), maggs_setup_local_lattice(), p3m_assign_charge(), and p3m_calc_local_ca_mesh().
| int node_grid[3] = { 0, 0, 0} |
The number of nodes in each spatial dimension.
Definition at line 48 of file grid.c.
Referenced by calc_node_neighbors(), calc_processor_min_num_cells(), dd_exchange_and_sort_particles(), dd_prepare_comm(), dd_update_communicators_w_boxl(), dfft_init(), dp3m_sanity_checks(), fft_init(), grid_changed_box_l(), grid_changed_n_nodes(), halo_push_communication(), init_lattice(), layered_topology_init(), lb_calc_densprof(), lb_calc_velprof(), maggs_calc_init_e_field(), maggs_sanity_checks(), map_lattice_to_node(), map_position_node_array(), mpi_init(), node_grid_is_set(), p3m_sanity_checks(), prepare_halo_communication(), tclcallback_node_grid(), tclcommand_analyze_fluid_parse_densprof(), and tclcommand_cellsystem().
| int node_neighbors[6] = {0, 0, 0, 0, 0, 0} |
the six nearest neighbors of a node in the node grid.
Definition at line 50 of file grid.c.
Referenced by calc_node_neighbors(), dd_exchange_and_sort_particles(), dd_prepare_comm(), dp3m_calc_send_mesh(), dp3m_gather_fft_grid(), dp3m_spread_force_grid(), grid_changed_n_nodes(), halo_push_communication(), init_lattice(), lb_check_halo_regions(), maggs_calc_init_e_field(), maggs_exchange_surface_patch(), p3m_calc_send_mesh(), p3m_gather_fft_grid(), p3m_spread_force_grid(), and prepare_halo_communication().
| int node_pos[3] = {-1,-1,-1} |
position of node in node grid
Definition at line 49 of file grid.c.
Referenced by calc_node_neighbors(), dd_exchange_and_sort_particles(), dd_prepare_comm(), dd_update_communicators_w_boxl(), dfft_init(), dp3m_calc_send_mesh(), dp3m_gather_fft_grid(), dp3m_spread_force_grid(), fft_init(), grid_changed_box_l(), grid_changed_n_nodes(), lb_calc_densprof(), lb_calc_velprof(), maggs_calc_init_e_field(), mpi_init(), p3m_calc_send_mesh(), p3m_gather_fft_grid(), and p3m_spread_force_grid().
| int periodic = 7 |
Flags for all three dimensions wether pbc are applied (default).
The first three bits give the periodicity
Definition at line 52 of file grid.c.
Referenced by mpi_init(), mpi_local_stress_tensor_slave(), tclcallback_periodicity(), and tclcommand_analyze_parse_local_stress_tensor().
1.7.5.1