Functions/Subroutines | |
program | diffusion_calc |
Program diffusion_calc. More... | |
subroutine | calculate_trans_msd (incoords, msd, max_steps, binned, bin_width, bin_offset, decomp, msd_decomp) |
Read in all coordinates and calculate the translational MSD. More... | |
subroutine | calculate_rot_msd (inorient, msd, max_steps, inheight, binned, bin_width, bin_offset) |
Read in all orientations and calculate the unbounded rotational MSD – See Mazza et al Phys. Rev. E (2007), 76, 031203. More... | |
subroutine | calculate_rot_acf (inorient, acf, max_steps, inheight, binned, bin_width, bin_offset) |
Read in all orientations and calculate the rotational autocorrelation function. More... | |
subroutine | write_function (io_data, funct, interval, binned, bin_width, bin_offset, desc) |
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 calculate_rot_acf | ( | real(kind=8), dimension(:,:,:,:), intent(in) | inorient, |
real(kind=8), dimension(:,:), intent(out) | acf, | ||
integer, intent(in) | max_steps, | ||
real(kind=8), dimension(:,:), intent(in) | inheight, | ||
logical | binned, | ||
real(kind=8), intent(in) | bin_width, | ||
real(kind=8), intent(in) | bin_offset | ||
) |
Read in all orientations and calculate the rotational autocorrelation function.
inorient | :: 4D array containing coordinates of each solute in each frame. Order (axis, axis vector, solute, frame) |
msd | :: 2D array for storing the ACF, allows height binning. Order (msd,bin) – MUST BE ALLOCATED IN CALLER |
max_steps | :: number of frames in the largest time interval over which ACF is calculated |
inheight | :: 2D vector giving the height of each solute in each frame. Needed for binning |
binned | :: is the ACF to be binned in z dimension |
bin_width | :: If binning data in z, what is the width in Ang of the bins |
bin_offset | :: If binning data in z, at what height does the first bin begin |
subroutine calculate_rot_msd | ( | real(kind=8), dimension(:,:,:,:), intent(in) | inorient, |
real(kind=8), dimension(:,:), intent(out) | msd, | ||
integer, intent(in) | max_steps, | ||
real(kind=8), dimension(:,:), intent(in) | inheight, | ||
logical | binned, | ||
real(kind=8), intent(in) | bin_width, | ||
real(kind=8), intent(in) | bin_offset | ||
) |
Read in all orientations and calculate the unbounded rotational MSD – See Mazza et al Phys. Rev. E (2007), 76, 031203.
Uses z axis vector as principle axis
inorient | :: 4D array containing coordinates of each solute in each frame. Order (axis, axis vector, solute, frame) |
msd | :: 2D array for storing the MSD, allows height binning. Order (msd,bin) – MUST BE ALLOCATED IN CALLER |
max_steps | :: number of frames in the largest time interval over which MSD is calculated |
inheight | :: 2D vector giving the height of each solute in each frame. Needed for binning |
binned | :: is the MSD to be binned in z dimension |
bin_width | :: If binning data in z, what is the width in Ang of the bins |
bin_offset | :: If binning data in z, at what height does the first bin begin |
subroutine calculate_trans_msd | ( | real(kind=8), dimension(:,:,:), intent(in) | incoords, |
real(kind=8), dimension(:,:), intent(out) | msd, | ||
integer, intent(in) | max_steps, | ||
logical | binned, | ||
real(kind=8), intent(in) | bin_width, | ||
real(kind=8), intent(in) | bin_offset, | ||
logical | decomp, | ||
real(kind=8), dimension(:,:,:), intent(out) | msd_decomp | ||
) |
Read in all coordinates and calculate the translational MSD.
incoords | :: 3D array containing coordinates of each solute in each frame. Order (x/y/z, solute, frame) |
msd | :: 2D array for storing the MSD, allows height binning. Order (msd,bin) – MUST BE ALLOCATED IN CALLER |
max_steps | :: number of frames in the largest time interval over which MSD is calculated |
binned | :: is the MSD to be binned in z dimension |
bin_width | :: If binning data in z, what is the width in Ang of the bins |
bin_offset | :: If binning data in z, at what height does the first bin begin |
decomp | :: Should diffusion be decomposed into 1D diffusion in z, and 2D diffusion in xy |
msd_decomp | :: 3D array for storing the decomposed MSD, allows height binning. Order ([xy,z],msd,bin) – MUST BE ALLOCATED IN CALLER |
subroutine write_function | ( | integer, intent(in) | io_data, |
real(kind=8), dimension(:,:), intent(in) | funct, | ||
real(kind=8), intent(in) | interval, | ||
logical, intent(in) | binned, | ||
real(kind=8), intent(in) | bin_width, | ||
real(kind=8), intent(in) | bin_offset, | ||
character(76), intent(in) | desc | ||
) |