ESPResSo 3.2.0-167-g2c9ead1-git
Extensible Simulation Package for Soft Matter Research
Defines
grid.h File Reference

Domain decomposition for parallel computing. More...

#include "utils.h"
#include <limits.h>
#include "communication.h"
#include "errorhandling.h"
Include dependency graph for grid.h:

Go to the source code of this file.

Defines

#define PERIODIC(coord)   (periodic & (1L << coord))
 Macro that tests for a coordinate being periodic or not.

Functions

Exported 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.
MDINLINE void map_node_array (int node, int pos[3])
 node mapping: array -> node.
MDINLINE int map_array_node (int pos[3])
 node mapping: node -> array.
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_n_nodes ()
 called from mpi_bcast_parameter .
void grid_changed_box_l ()
 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
MDINLINE void get_mi_vector (double res[3], double a[3], double b[3])
 get the minimal distance vector of two vectors in the current bc.
MDINLINE void fold_coordinate (double pos[3], int image_box[3], int dir)
 fold a coordinate to primary simulation box.
MDINLINE void fold_position (double pos[3], int image_box[3])
 fold particle coordinates to primary simulation box.
MDINLINE void unfold_position (double pos[3], int image_box[3])
 unfold coordinates to physical position.

Variables

Exported Variables
int node_grid [3]
 The number of nodes in each spatial dimension.
int node_pos [3]
 position of node in node grid
int node_neighbors [6]
 the six nearest neighbors of a node in the node grid.
int boundary [6]
 where to fold particles that leave local box in direction i.
int periodic
 Flags for all three dimensions wether pbc are applied (default).
double box_l [3]
 Simulation box dimensions.
double box_l_i [3]
 1 / box dimensions.
double min_box_l
 Smallest simulation box dimension (box_l).
double local_box_l [3]
 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]
 Left (bottom, front) corner of this nodes local box.
double my_right [3]
 Right (top, back) corner of this nodes local box.

Detailed Description

Domain decomposition for parallel computing.

The primary simulation box is divided into orthogonal rectangular subboxes which are assigned to the different nodes (or processes or threads if you want). This grid is described in node_grid. Each node has a number this_node and a position node_pos in that grid. Each node has also 6 nearest neighbors node_neighbors which are necessary for the communication between the nodes (see also ghosts.c and p3m.c for more details about the communication.

For the 6 directions we have the following convention:

directions.gif
Convention for the order of the directions

The Figure illustrates the direction convetion used for arrays with 6 (e.g. node_neighbors, boundary) and 3 entries (e.g node_grid, box_l , my_left,...).

For more information on the domain decomposition, see grid.c.

Definition in file grid.h.


Define Documentation

#define PERIODIC (   coord)    (periodic & (1L << coord))

Function Documentation

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.

Parameters:
nodenumber 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().

MDINLINE void fold_coordinate ( double  pos[3],
int  image_box[3],
int  dir 
)

fold a coordinate to primary simulation box.

Parameters:
posthe position...
image_boxand the box
dirthe coordinate to fold: dir = 0,1,2 for x, y and z coordinate.

Both pos and image_box are I/O, i. e. a previously folded position will be folded correctly.

Definition at line 197 of file grid.h.

References box_l, box_l_i, ERROR_SPRINTF, ES_DOUBLE_SPACE, ES_INTEGER_SPACE, PERIODIC, ROUND_ERROR_PREC, and runtime_error().

Referenced by calc_diffusion_profile(), calc_radial_density_map(), dd_append_particles(), dd_exchange_and_sort_particles(), density_profile_av(), fold_all(), and fold_position().

MDINLINE void fold_position ( double  pos[3],
int  image_box[3] 
)
MDINLINE void get_mi_vector ( double  res[3],
double  a[3],
double  b[3] 
)

get the minimal distance vector of two vectors in the current bc.

Parameters:
athe vector to subtract from
bthe vector to subtract
reswhere to store the result

Definition at line 176 of file grid.h.

References box_l, box_l_i, dround(), i, and PERIODIC.

Referenced by add_bonded_energy(), add_bonded_force(), add_bonded_virials(), add_ljangle_pair_force(), add_three_body_bonded_stress(), angle_cosine_energy(), angle_cossquare_energy(), angle_energy(), angle_harmonic_energy(), angledist_energy(), bilayer_density_profile_sphere(), calc_angle_3body_forces(), calc_angle_3body_tabulated_forces(), calc_angle_cosine_3body_forces(), calc_angle_cosine_force(), calc_angle_cossquare_3body_forces(), calc_angle_cossquare_force(), calc_angle_force(), calc_angle_harmonic_3body_forces(), calc_angle_harmonic_force(), calc_angledist_force(), calc_angledist_param(), calc_area_force_local_complicated(), calc_dihedral_angle(), calc_dipole_dipole_ia(), calc_dipole_of_molecule(), calc_endangledist_pair_force(), calc_force_between_mol(), calc_mol_pos(), calc_overlap_angle_force(), calc_pwangle(), calc_scaling(), calc_scaling2(), calc_tab_angle_force(), calc_wallbondyz(), calc_wallrdfyz(), compute_pos_corr_vec(), compute_vel_corr_vec(), distribute_mol_force(), distribute_tensors(), distto(), does_line_go_through_cube(), ELC_P3M_dielectric_layers_energy_contribution(), ELC_P3M_dielectric_layers_energy_self(), ELC_P3M_dielectric_layers_force_contribution(), ELC_P3M_self_forces(), get_mol_dist(), get_mol_dist_vector_from_molid_cfg(), get_nonbonded_interaction(), incubewithskin(), integrate_reaction(), ljangle_pair_energy(), local_stress_tensor_calc(), min_distance2(), nbhood(), nsq_calculate_energies(), nsq_calculate_ia(), nsq_calculate_ia_iccp3m(), nsq_calculate_virials(), observable_interacts_with(), overlap_angle_energy(), print_bond_len(), reducepos(), tab_angle_energy(), test_mesh_elements(), and vs_relate_to().

void grid_changed_box_l ( )
void grid_changed_n_nodes ( )
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.

Parameters:
g3d3d grid.
g2d2d grid.
multfactors between 3d and 2d grid dimensions
Returns:
index of the row direction [0,1,2].

Definition at line 196 of file grid.c.

References i.

Referenced by dfft_init(), and fft_init().

MDINLINE int map_array_node ( int  pos[3])

node mapping: node -> array.

Returns:
rank of the node at position pos.
Parameters:
posposition of the node in node grid.

Definition at line 120 of file grid.h.

References comm_cart, and MPI_Cart_rank().

Referenced by map_lattice_to_node(), and map_position_node_array().

MDINLINE void map_node_array ( int  node,
int  pos[3] 
)

node mapping: array -> node.

Parameters:
noderank of the node you want to know the position for.
posposition of the node in node grid.

Definition at line 110 of file grid.h.

References comm_cart, and MPI_Cart_coords().

Referenced by calc_node_neighbors(), dfft_init(), fft_init(), lb_calc_densprof(), lb_calc_velprof(), and lb_init_boundaries().

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

return wether node grid was set.

Definition at line 70 of file grid.c.

References node_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

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

MDINLINE void unfold_position ( double  pos[3],
int  image_box[3] 
)

unfold coordinates to physical position.

Parameters:
posthe position...
image_boxand the box

Both pos and image_box are I/O, i.e. image_box will be (0,0,0) afterwards.

Definition at line 237 of file grid.h.

References box_l, and i.

Referenced by calc_local_mol_info(), local_system_CMS(), meta_perform(), tclcommand_analyze_parse_and_print_check_mol(), tclcommand_part_print_position(), tclcommand_writemd(), and updatePartCfg().


Variable Documentation

int boundary[6]
double box_l[3]

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(), on_coulomb_change(), 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]
double local_box_l[3]
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().

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]
double my_right[3]
int node_grid[3]
int node_pos[3]
int periodic

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