Compute the radial distribution function. More...
Functions/Subroutines | |
program | rdf |
Program rdf. More... | |
subroutine | calc_frame_rdf (rdf_hist, position, start_centre, end_centre, start_dist, end_dist, box_dim) |
Calculate the rdf of a single trajectory frame. More... | |
subroutine | normalise_rdf (unnormalised_rdf, normalised_rdf, counter, nb_centre, dens) |
Normalise radial distribution function. More... | |
subroutine | virial (rdf_hist, virial_integral, rdf_integral) |
Subroutine find osmotic second virial coefficients by integration of RDF. More... | |
subroutine | write_rdf_log (log_fileno, rdf_hist, virial_int, rdf_int) |
Subroutine to write normalised RDF to file. More... | |
Compute the radial distribution function.
Copyright (c) 2009, 2010, 2015, 2016, 2019 Heidelberg Institute of Theoretical Studies (HITS, www.h-its.org) Schloss-Wolfsbrunnenweg 35 69118 Heidelberg, Germany
Please send your contact address to get information on updates and new features to "mcmsoft@h-its.org". Questions will be answered as soon as possible.
References: see also http://mcm.h-its.org/sda7/do:c/doc_sda7/references.html:
Brownian dynamics simulation of protein-protein diffusional encounter. (1998) Methods, 14, 329-341.
SDA 7: A modular and parallel implementation of the simulation of diffusional association software. Journal of computational chemistry 36.21 (2015): 1631-1645.
Authors: M.Martinez, N.J.Bruce, J.Romanowska, D.B.Kokh, P.Mereghetti, X. Yu, M. Ozboyaci, M. Reinhardt, P. Friedrich, R.R.Gabdoulline, S.Richter and R.C.Wade
subroutine calc_frame_rdf | ( | type(type_int_hist), intent(inout) | rdf_hist, |
real(kind=8), dimension(:,:), intent(in) | position, | ||
integer, intent(in) | start_centre, | ||
integer, intent(in) | end_centre, | ||
integer, intent(in) | start_dist, | ||
integer, intent(in) | end_dist, | ||
real(kind=8), dimension(3), intent(in) | box_dim | ||
) |
Calculate the rdf of a single trajectory frame.
rdf_hist | : instance of mod_histogram::type_int_hist used to store RDF |
tab_prot | : array to store the positions of all solutes in current frame |
start_centre | : Index to first solute of central type |
end_centre | : Index to last solute of central type |
start_dist | : Index to first solute of distribution type |
end_dist | : Index to last solute of distribution type |
box_dim | : Array to store xyz box dimensions |
subroutine normalise_rdf | ( | type(type_int_hist), intent(in) | unnormalised_rdf, |
type(type_real_hist), intent(out) | normalised_rdf, | ||
integer, intent(in) | counter, | ||
integer, intent(in) | nb_centre, | ||
real(kind=8), intent(in) | dens | ||
) |
Normalise radial distribution function.
unnormalised_rdf | : mod_histogram::type_int_hist instance used to store input RDF |
normalised_rdf | : mod_histogram::type_real_hist instance used to store normalised output RDF |
counter | : Number of frames used in generating distribution |
nb_centre | : Number of central solutes used to create distribution |
protein_array | : Needed to calculate number density of each solute |
subroutine virial | ( | type(type_real_hist), intent(in) | rdf_hist, |
real(kind=8), intent(out) | virial_integral, | ||
real(kind=8), intent(out) | rdf_integral | ||
) |
Subroutine find osmotic second virial coefficients by integration of RDF.
rdf_hist | : normalised real RDF histogram |
virial_integral | : integral \( \int \left( g\left(r\right) - 1\right)r^2 \textrm{d}r \) over histogram length |
rdf_integral | : integral \( \int \left( g\left(r\right) - 1\right) \textrm{d}r \) over histogram length |
Integrates bins by midpoint rectangular quadrature
subroutine write_rdf_log | ( | integer, intent(in) | log_fileno, |
type(type_real_hist), intent(in) | rdf_hist, | ||
real(kind=8), intent(in) | virial_int, | ||
real(kind=8), intent(in) | rdf_int | ||
) |
Subroutine to write normalised RDF to file.
log_fileno | : File number to write log file. File already opened in calling routine |
rdf_hist | : Normalised RDF stored in mod_histogram::type_real_hist instance |