(git:6a2e663)
mc_coordinates Module Reference

contains miscellaneous subroutines used in the Monte Carlo runs,mostly geared towards changes in coordinates More...

Functions/Subroutines

subroutine, public check_for_overlap (force_env, nchains, nunits, loverlap, mol_type, cell_length, molecule_number)
 looks for overlaps (intermolecular distances less than rmin) More...
 
subroutine, public get_center_of_mass (coordinates, natom, center_of_mass, mass)
 calculates the center of mass of a given molecule More...
 
subroutine, public mc_coordinate_fold (coordinates, nchains_tot, mol_type, mass, nunits, box_length)
 folds all the coordinates into the center simulation box using a center of mass cutoff More...
 
subroutine, public generate_cbmc_swap_config (force_env, BETA, max_val, min_val, exp_max_val, exp_min_val, nswapmoves, rosenbluth_weight, start_atom, natoms_tot, nunits, nunits_mol, mass, loverlap, choosen_energy, old_energy, ionode, lremove, mol_type, nchains, source, group, rng_stream, avbmc_atom, rmin, rmax, move_type, target_atom)
 takes the last molecule in a force environment and moves it around to different center of mass positions and orientations, selecting one based on the rosenbluth weight More...
 
subroutine, public rotate_molecule (r, mass, natoms, rng_stream)
 rotates a molecule randomly around the center of mass, sequentially in x, y, and z directions More...
 
subroutine, public find_mc_test_molecule (mc_molecule_info, start_atom, box_number, molecule_type, rng_stream, box, molecule_type_old)
 selects a molecule at random to perform a MC move on...you can specify the box the molecule should be in, its type, both, or neither More...
 
subroutine, public create_discrete_array (cell, discrete_array, step_size)
 generates an array that tells us which sides of the simulation cell we can increase or decrease using a discrete volume move More...
 
subroutine, public cluster_search (mc_par, force_env, cluster, nchains, nunits, mol_type, total_clus)
 determine the number of cluster present in the given configuration based on the rclus value More...
 

Detailed Description

contains miscellaneous subroutines used in the Monte Carlo runs,mostly geared towards changes in coordinates

Author
MJM

Function/Subroutine Documentation

◆ check_for_overlap()

subroutine, public mc_coordinates::check_for_overlap ( type(force_env_type), pointer  force_env,
integer, dimension(:), intent(in)  nchains,
integer, dimension(:), intent(in)  nunits,
logical, intent(out)  loverlap,
integer, dimension(:), intent(in)  mol_type,
real(kind=dp), dimension(1:3), intent(in), optional  cell_length,
integer, intent(in), optional  molecule_number 
)

looks for overlaps (intermolecular distances less than rmin)

Parameters
force_envthe force environment containing the coordinates
nchainsthe number of molecules of each type in the box
nunitsthe number of interaction sites for each molecule
loverlapreturns .TRUE. if atoms overlap
mol_typean array that indicates the type of each molecule
cell_lengththe length of the box...if none is specified, it uses the cell found in the force_env
molecule_numberif present, just look for overlaps with this molecule

Suitable for parallel use.

Author
MJM

Definition at line 67 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_center_of_mass()

subroutine, public mc_coordinates::get_center_of_mass ( real(kind=dp), dimension(:, :), intent(in)  coordinates,
integer, intent(in)  natom,
real(kind=dp), dimension(1:3), intent(out)  center_of_mass,
real(kind=dp), dimension(:), intent(in)  mass 
)

calculates the center of mass of a given molecule

Parameters
coordinatesthe coordinates of the atoms in the molecule
natomthe number of atoms in the molecule
center_of_massthe coordinates of the center of mass
massthe mass of the atoms in the molecule

Designed for parallel use.

Author
MJM

Definition at line 188 of file mc_coordinates.F.

Here is the caller graph for this function:

◆ mc_coordinate_fold()

subroutine, public mc_coordinates::mc_coordinate_fold ( real(kind=dp), dimension(:, :), intent(inout)  coordinates,
integer, intent(in)  nchains_tot,
integer, dimension(:), intent(in)  mol_type,
real(kind=dp), dimension(:, :), intent(in)  mass,
integer, dimension(:), intent(in)  nunits,
real(kind=dp), dimension(1:3), intent(in)  box_length 
)

folds all the coordinates into the center simulation box using a center of mass cutoff

Parameters
coordinatesthe coordinates of the atoms in the system
nchains_totthe total number of molecules in the box
mol_typean array that indicates the type of every molecule in the box
massthe mass of every atom for all molecule types
nunitsthe number of interaction sites for each molecule type
box_lengthan array for the lengths of the simulation box sides
 Designed for parallel use.
Author
MJM

Definition at line 235 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ generate_cbmc_swap_config()

subroutine, public mc_coordinates::generate_cbmc_swap_config ( type(force_env_type), pointer  force_env,
real(kind=dp), intent(in)  BETA,
real(kind=dp), intent(in)  max_val,
real(kind=dp), intent(in)  min_val,
real(kind=dp), intent(in)  exp_max_val,
real(kind=dp), intent(in)  exp_min_val,
integer, intent(in)  nswapmoves,
real(kind=dp), intent(out)  rosenbluth_weight,
integer, intent(in)  start_atom,
integer, intent(in)  natoms_tot,
integer, dimension(:), intent(in)  nunits,
integer, intent(in)  nunits_mol,
real(dp), dimension(1:nunits_mol), intent(in)  mass,
logical, intent(out)  loverlap,
real(kind=dp), intent(out)  choosen_energy,
real(kind=dp), intent(in)  old_energy,
logical, intent(in)  ionode,
logical, intent(in)  lremove,
integer, dimension(:), intent(in)  mol_type,
integer, dimension(:), intent(in)  nchains,
integer, intent(in)  source,
type(mp_comm_type)  group,
type(rng_stream_type), intent(inout)  rng_stream,
integer, intent(in), optional  avbmc_atom,
real(kind=dp), intent(in), optional  rmin,
real(kind=dp), intent(in), optional  rmax,
character(len=*), intent(in), optional  move_type,
integer, intent(in), optional  target_atom 
)

takes the last molecule in a force environment and moves it around to different center of mass positions and orientations, selecting one based on the rosenbluth weight

Parameters
force_envthe force environment containing the coordinates
BETAthe value of 1/kT for this simulations, in a.u.
max_val...
min_val...
exp_max_val...
exp_min_val...
nswapmovesthe number of desired trial configurations
rosenbluth_weightthe Rosenbluth weight for this set of configs
start_atomthe atom number that the molecule to be swapped starts on
natoms_totthe total number of interaction sites in the box
nunitsthe number of interaction sites for every molecule_type
nunits_mol...
massthe mass for every interaction site of every molecule type
loverlapthe flag that indicates if all of the configs have an atomic overlap
choosen_energythe energy of the chosen config
old_energythe energy that we subtract from all of the trial energies to prevent numerical overflows
ionodeindicates if we're on the main CPU
lremoveis this the Rosenbluth weight for a removal box?
mol_typean array that contains the molecule type for every atom in the box
nchainsthe number of molecules of each type in this box
sourcethe MPI source value, for broadcasts
groupthe MPI group value, for broadcasts
rng_streamthe random number stream that we draw from
avbmc_atom...
rmin...
rmax...
move_type...
target_atom...
Optional Avbmc Flags
  • avbmc_atom: the atom number that serves for the target atom in each molecule (1 is the first atom in the molecule, etc.)
  • rmin: the minimum AVBMC radius for the shell around the target
  • rmax: the maximum AVBMC radius for the shell around the target
  • move_type: generate configs in the "in" or "out" volume
  • target_atom: the number of the avbmc atom in the target molecule
Suitable for parallel.
Author
MJM

Definition at line 326 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rotate_molecule()

subroutine, public mc_coordinates::rotate_molecule ( real(kind=dp), dimension(1:3, 1:natoms), intent(inout)  r,
real(kind=dp), dimension(1:natoms), intent(in)  mass,
integer, intent(in)  natoms,
type(rng_stream_type), intent(inout)  rng_stream 
)

rotates a molecule randomly around the center of mass, sequentially in x, y, and z directions

Parameters
rthe coordinates of the molecule to rotate
massthe mass of all the atoms in the molecule
natomsthe number of atoms in the molecule
rng_streamthe stream we pull random numbers from

Use only in serial.

Author
MJM

Definition at line 626 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ find_mc_test_molecule()

subroutine, public mc_coordinates::find_mc_test_molecule ( type(mc_molecule_info_type), pointer  mc_molecule_info,
integer, intent(out)  start_atom,
integer, intent(out)  box_number,
integer, intent(out)  molecule_type,
type(rng_stream_type), intent(inout)  rng_stream,
integer, intent(in), optional  box,
integer, intent(in), optional  molecule_type_old 
)

selects a molecule at random to perform a MC move on...you can specify the box the molecule should be in, its type, both, or neither

Parameters
mc_molecule_infothe structure that contains some global molecule data
start_atomthe number of the first atom in the chosen molecule in relation to the force_env it's in
box_numberthe box the chosen molecule is in
molecule_typethe type of molecule the chosen molecule is
rng_streamthe stream we pull random numbers from
boxif present, tells the routine which box to grab a molecule from
molecule_type_oldif present, tells the routine which molecule type to select from
Author
MJM

Definition at line 714 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ create_discrete_array()

subroutine, public mc_coordinates::create_discrete_array ( real(dp), dimension(1:3), intent(in)  cell,
integer, dimension(1:3, 1:2), intent(out)  discrete_array,
real(dp), intent(in)  step_size 
)

generates an array that tells us which sides of the simulation cell we can increase or decrease using a discrete volume move

Parameters
cellthe lengths of the sides of the cell
discrete_arraythe array that indicates which sides we can move
step_sizethe size of the discrete volume move

Suitable for parallel.

Author
MJM

Definition at line 866 of file mc_coordinates.F.

Here is the caller graph for this function:

◆ cluster_search()

subroutine, public mc_coordinates::cluster_search ( type(mc_simpar_type), pointer  mc_par,
type(force_env_type), pointer  force_env,
integer, dimension(:, :), intent(inout)  cluster,
integer, dimension(:), intent(in)  nchains,
integer, dimension(:), intent(in)  nunits,
integer, dimension(:), intent(in)  mol_type,
integer, intent(inout)  total_clus 
)

determine the number of cluster present in the given configuration based on the rclus value

Parameters
mc_parthe mc parameters for the force env
force_envthe force environment containing the coordinates
cluster...
nchainsthe number of molecules of each type in the box
nunitsthe number of interaction sites for each molecule
mol_typean array that indicates the type of each molecule
total_clus...
Original Multiparticle/Cluster Translation paper: Orkoulas, Gerassimos, and Athanassios Z. Panagiotopoulos. Free energy and phase equilibria for the restricted primitive model of ionic fluids from Monte Carlo simulations. J. Chem. Phys. 1994,101.2,1452-1459.
Author
Himanshu Goel

Definition at line 1013 of file mc_coordinates.F.

Here is the call graph for this function:
Here is the caller graph for this function: