| Functions/Subroutines | |
| subroutine | print_acf (ACF, time_snapshot, b_type, value_type, title) | 
| General function for printing ACFs type, in a single file.  More... | |
| subroutine | print_coeffd_msd (D_COEFF, region_fit, b_type, value_type, dt, b_trans, b_deriv, io_msd) | 
| Print coefficient from the msd fit method.  More... | |
| subroutine | fit_msd_by_regions (input, nb_region, dt, b_trans, io_dt) | 
| Fit the MSD graph.  More... | |
| subroutine | fit_msd_derivative (input, D_COEFF, shift_step, array_limit) | 
| Fit the derivative of the MSD.  More... | |
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
Common functions used by diffusion program tools
| subroutine fit_msd_by_regions | ( | real ( kind=8 ), dimension(:), intent(in) | input, | 
| integer, intent(in) | nb_region, | ||
| real ( kind=4 ), intent(in) | dt, | ||
| logical, intent(in) | b_trans, | ||
| integer, intent(in) | io_dt | ||
| ) | 
Fit the MSD graph.
| input | : 1 dimension array | 
| nb_region | : how many blocks are divided | 
| dt | : time between each frame | 
| io_dt | : file unit of the output file. Should be opened / close by calller | 
| b_trans | : parameter trans / rot for factor 4 or 6 | 
| subroutine fit_msd_derivative | ( | real ( kind=8 ), dimension(:), intent(in) | input, | 
| real ( kind=8 ), dimension(4) | D_COEFF, | ||
| integer | shift_step, | ||
| integer, dimension(4) | array_limit | ||
| ) | 
Fit the derivative of the MSD.
Different method
 To add array for storing limit of fit, will be easier to print out 
| input | : 1D array with values to fit | 
| D_COEFF | : array to store the results | 
| shift_step | : ?? | 
| array_limit | : store the fitted region | 
| subroutine print_acf | ( | real ( kind=8 ), dimension (:,:,:), intent(in) | ACF, | 
| real ( kind=4 ), intent(in) | time_snapshot, | ||
| logical, intent(in) | b_type, | ||
| integer, intent(in) | value_type, | ||
| character, dimension(*), intent(in) | title | ||
| ) | 
General function for printing ACFs type, in a single file.
| ACF | : 3 diemsional array for autocorrelation function ( "x,y,z", protein_nb, time ) | 
| time_snapshot | : timestep between the frames | 
| b_type | : if true, there is different types of solutes in the simulation | 
| value_type | : indicates which type if b_type is true | 
| title | : filename for the output | 
| subroutine print_coeffd_msd | ( | real ( kind=8 ), dimension (:,:,:) | D_COEFF, | 
| integer, dimension(:,:) | region_fit, | ||
| logical | b_type, | ||
| integer | value_type, | ||
| real ( kind=4 ) | dt, | ||
| logical | b_trans, | ||
| logical | b_deriv, | ||
| integer | io_msd | ||
| ) | 
Print coefficient from the msd fit method.
Give io_msd to append the data to the main log file
| D_COEFF | : contains the computed values | 
| region_fit | : info about the fitting | 
| b_type | : if true, there is different types of solutes in the simulation | 
| value_type | : indicates which type if b_type is true | 
| dt | : timestep | 
| b_trans | : true if translation, false rotation | 
| b_deriv | : true if exact derivation, false simplified derivative ( MSD/time ) | 
| io_msd | : file descriptor | 
 1.8.13
Imprint/Privacy
 1.8.13
Imprint/Privacy