SDA (SDA flex)  7.2
Simulation of Diffusional Association
Data Types | Functions/Subroutines
mod_rate_calculation Module Reference

Define type association_rate. More...

Data Types

type  association_rate
 General type used storing statistics about the runs. More...
 

Functions/Subroutines

subroutine deallocate_rate_association (this)
 Delete arrays and reinitialize data members. More...
 
subroutine update_escape (this, nb_traj, visited_windows, escape, sum_time, first_time)
 Update the arrays after each trajectory. More...
 
subroutine allocate_rate_association (this, first_window, nb_window, width_window, nb_contact, bsurf, csurf, dm, ibootsp, ifpt, iboot_et, stop_traj, anal_corr, dh_rad_sum, real_net_charge_product, ionic)
 Initialize association_rate and allocate arrays. More...
 
subroutine calc_rate_corr (brate, crate, bsurf, csurf, dh_rad_sum, real_net_charge_product, ionic, dm)
 
subroutine write_boot_file (this, total_run)
 Write the bootstrap output files. More...
 
subroutine print_table_rate (this, total_run, opt_koff)
 Print final results, could be used for restart or checkpoint. More...
 
subroutine print_sumtime (this, norm)
 Print only the array_sumtime, implemented for sdamm and association_rate. More...
 

Detailed Description

Define type association_rate.

Used only with sda_2proteins

Function/Subroutine Documentation

◆ allocate_rate_association()

subroutine mod_rate_calculation::allocate_rate_association ( type ( association_rate this,
real ( kind=4 ), intent(in)  first_window,
integer, intent(in)  nb_window,
real ( kind=4 ), intent(in)  width_window,
integer, intent(in)  nb_contact,
real ( kind=4 ), intent(in)  bsurf,
real ( kind=4 ), intent(in)  csurf,
real ( kind=8 ), intent(in)  dm,
logical  ibootsp,
logical  ifpt,
logical  iboot_et,
logical  stop_traj,
logical  anal_corr,
real ( kind=8 ), intent(in)  dh_rad_sum,
real ( kind=8 ), intent(in)  real_net_charge_product,
real ( kind=4 ), intent(in)  ionic 
)

Initialize association_rate and allocate arrays.

dm need to be the sum of the diffusion coeff.

Parameters
this: instance of association_rate
first_window: distance of the first window, in Angstrom
nb_window: number of windows
width_window: width of each window
nb_contact: number of contacts to record, from 1 to nb_contact included
bsurf,csurf: values of the b and c surface
dm: use the sum of the diffusion terms, absolute trans. diff. coeff.
ibootsp,ifpt,iboot_et: file descriptor for the bootstrap output files
stop_traj: should trajectory end when most strict reaction criteria met
anal_corr: should analytical correction be applied to account for interactions at b surface
dh_rad_sum: sum of Debye-Hueckel cavity radii of two solutes
real_net_charge_product: product of net charge of two solutes
ionic: ionic strength
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_rate_corr()

subroutine mod_rate_calculation::calc_rate_corr ( real(kind=8), intent(out)  brate,
real(kind=8), intent(out)  crate,
real(kind=4), intent(in)  bsurf,
real(kind=4), intent(in)  csurf,
real(kind=8), intent(in)  dh_rad_sum,
real(kind=8), intent(in)  real_net_charge_product,
real(kind=4), intent(in)  ionic,
real(kind=8), intent(in)  dm 
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deallocate_rate_association()

subroutine mod_rate_calculation::deallocate_rate_association ( type ( association_rate this)

Delete arrays and reinitialize data members.

Parameters
this: instance of association_rate
Here is the caller graph for this function:

◆ print_sumtime()

subroutine mod_rate_calculation::print_sumtime ( type ( association_rate this,
real ( kind=4 )  norm 
)

Print only the array_sumtime, implemented for sdamm and association_rate.

Parameters
this: instance of association_rate
norm: normalisation factor
Here is the call graph for this function:

◆ print_table_rate()

subroutine mod_rate_calculation::print_table_rate ( type ( association_rate this,
integer  total_run,
logical, optional  opt_koff 
)

Print final results, could be used for restart or checkpoint.

Print data in the main output file
opt_koff, will not print the rate, does not make sense
if sda_2proteins, print all

Parameters
this: instance of association_rate
total_run: total number of runs
opt_koff: optionaly print association rates

◆ update_escape()

subroutine mod_rate_calculation::update_escape ( type ( association_rate this,
integer, intent(in)  nb_traj,
integer, dimension(:), intent(in)  visited_windows,
real ( kind=8 ), dimension(:,:), intent(in)  escape,
real ( kind=4 ), dimension (:,:), intent(in)  sum_time,
real ( kind=4 ), dimension (:,:), intent(in)  first_time 
)

Update the arrays after each trajectory.

Update all arrays at the end of trajectories
Function is thread-safe, with the use of ATOMIC sections minimum overhead.

Parameters
this: instance of association_rate
nb_traj: the total number of trajectories
visited_windows: array for association rate ( version with 4 integers )
escape: array for electron transfert, cannot anymore share the format with association rate
sum_time,first_time: external arrays belonging to each thread for updating global sumtime and firstime

◆ write_boot_file()

subroutine mod_rate_calculation::write_boot_file ( type ( association_rate this,
integer  total_run 
)

Write the bootstrap output files.

Call only at the end of the computation,

Parameters
this: instance of association_rate
total_run: the total number of runs really done (case of intercepted signal)
Imprint/Privacy