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

Module to define reaction criteria
More informations on the reaction criteria page. More...

Data Types

type  react_criter
 react_criter deals mainly with the reaction criteria positions and distances More...
 

Functions/Subroutines

subroutine allocate_react_criter (this, type_reaction_string, filename, et1f, et2f, nnons, nwrec, sdamd)
 Allocate arrays. More...
 
subroutine delete_react_criter (this)
 Delete the arrays. More...
 
subroutine allocate_rate (this, win0, nwin, dwin, nb_contact, b, c, dm, ibootsp, ifpt, stop_traj, anal_corr, stoke_rad_sum, real_net_charge_product, ionic)
 Initilialize the association_rate pointer. More...
 
subroutine read_rxna (this, filename)
 Read the reaction criteria file, usually *.rxna file. More...
 
subroutine read_et_rxna (this, et1f, et2f)
 Read reaction criteria for electron transfert calculation. More...
 
subroutine make_independent (this, dind)
 Make the list of independent pairs. More...
 
subroutine check_add_criteria (this, counter, nb_atom, xv, tmp_rcriter_pos, status_incl)
 Test when reading reaction criteria, called by read_pdb. More...
 
subroutine write_to_logfile (this, io_unit)
 Write information to a log file. More...
 

Variables

integer, parameter docking = 1
 Enumeration, associated with type_rection. More...
 
integer, parameter assoc_rate = 2
 
integer, parameter elec_transfer = 3
 
integer, parameter all = 4
 
integer, parameter off_rc = 5
 

Detailed Description

Module to define reaction criteria
More informations on the reaction criteria page.

Function/Subroutine Documentation

◆ allocate_rate()

subroutine mod_reactioncriteria::allocate_rate ( type ( react_criter this,
real ( kind=4 ), intent(in)  win0,
integer, intent(in)  nwin,
real ( kind=4 ), intent(in)  dwin,
integer, intent(in)  nb_contact,
real ( kind=4 ), intent(in)  b,
real ( kind=4 ), intent(in)  c,
real ( kind=8 ), intent(in)  dm,
logical  ibootsp,
logical  ifpt,
logical  stop_traj,
logical  anal_corr,
real ( kind=8 ), intent(in), optional  stoke_rad_sum,
real ( kind=8 ), intent(in), optional  real_net_charge_product,
real ( kind=4 ), intent(in), optional  ionic 
)

Initilialize the association_rate pointer.

Create the object and initiliaze the internal pointer p_assoc_rate

Parameters
this: instance of react_criter
win0,nwin,dwin: first_window, number and width of windows
nb_contact: the maximum number of contact to record
b,c: b and c surface
dm: translational diffusion coefficient of the solute 2
ibootsp,ifpt: if bootstrap for association rate or first tpassage time must be applied
stop_traj: should current traj stop when most strict reaction criteria met
Here is the call graph for this function:
Here is the caller graph for this function:

◆ allocate_react_criter()

subroutine mod_reactioncriteria::allocate_react_criter ( type ( react_criter this,
character*128, intent(in)  type_reaction_string,
character*128, intent(in)  filename,
character*128, intent(in)  et1f,
character*128, intent(in)  et2f,
integer, intent(in)  nnons,
integer, intent(in)  nwrec,
integer, intent(in)  sdamd 
)

Allocate arrays.

It sets up parameters and calls the function read_rxna or read_et_rxna for reading the reaction criteria The association_rate pointer may be then initilized with allocate_rate

Parameters
this: instance of react_criter
type_reaction_string: "docking" / "asociation" / "electron_transfert" / "all" / "off_rc"
filename: name of the file with the reaction criteria
et1f,et2f: in case of electron transfert the reaction criteria are defined in theses files
nnons: minimum number of nonspecific criteria to valid for accepting complexes
nwrec: in case of association rate, the specific window where reaction criteria are checked
Todo:
check initialization of dind
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_add_criteria()

subroutine mod_reactioncriteria::check_add_criteria ( type ( react_criter this,
integer  counter,
integer, intent(in)  nb_atom,
real ( kind=8 ), dimension( 3 ), intent(in)  xv,
real ( kind=8 ), dimension(:,:), pointer  tmp_rcriter_pos,
logical, intent(out)  status_incl 
)

Test when reading reaction criteria, called by read_pdb.

Used when one solute is flexible, search for the correct atom number
Check than the atom number of this particular conformation corresponds to the atom number of the initial reaction criteria

Parameters
this: instance of react_criter
counter: identifier of the solute
nb_atom: atom number
xv: position of the atom
tmp_rcriter_pos: return the position of the atom corresponding to the atom number
status_incl: return true if the atom has been included in the array of reaction criteria
Here is the caller graph for this function:

◆ delete_react_criter()

subroutine mod_reactioncriteria::delete_react_criter ( type ( react_criter this)

Delete the arrays.

Delete all allocated arrays and the association_rate object

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

◆ make_independent()

subroutine mod_reactioncriteria::make_independent ( type ( react_criter this,
real (kind=4), intent(in)  dind 
)

Make the list of independent pairs.

Todo:
Need to think about independent pairs in case of flexible/multi-proteins

Independent pairs apply only to associaiton rates
independence is check only for each pos_atom 1 and 2 separately but stored in matrix ( 1, 2 ) => no problem of center of mass

Parameters
this: instance of react_criter
dind: distance to use for indepence ( in Angstrom ), Default 6 A
Here is the caller graph for this function:

◆ read_et_rxna()

subroutine mod_reactioncriteria::read_et_rxna ( type ( react_criter this,
character*128, intent(in)  et1f,
character*128, intent(in)  et2f 
)

Read reaction criteria for electron transfert calculation.

Called by allocate_react_criter

Parameters
this: instance of react_criter
et1f,et2f: file names where the reaction criteria are defined
Here is the caller graph for this function:

◆ read_rxna()

subroutine mod_reactioncriteria::read_rxna ( type ( react_criter this,
character*128, intent(in)  filename 
)

Read the reaction criteria file, usually *.rxna file.

Called by allocate_react_criter
Do not check for independence, default ind_pairs is false , done optionaly by make_independent
Should be called the first time, if nb criteria are identical use update rxna
Data stored in arrays drc and drc2 in the order 1,nspec and nspec+1, nspec+1+aspec,..nspec+aspec+1+nnonspec

Parameters
this: instance of react_criter
filename: name of the reaction criteria file
Here is the caller graph for this function:

◆ write_to_logfile()

subroutine mod_reactioncriteria::write_to_logfile ( type ( react_criter this,
integer, intent(in)  io_unit 
)

Write information to a log file.

io_unit is given by read_input_file, egal to fort.52

Parameters
this: instance of react_criter
io_unit: file descriptor of the log file
Here is the caller graph for this function:

Variable Documentation

◆ all

integer, parameter mod_reactioncriteria::all = 4

◆ assoc_rate

integer, parameter mod_reactioncriteria::assoc_rate = 2

◆ docking

integer, parameter mod_reactioncriteria::docking = 1

Enumeration, associated with type_rection.

◆ elec_transfer

integer, parameter mod_reactioncriteria::elec_transfer = 3

◆ off_rc

integer, parameter mod_reactioncriteria::off_rc = 5
Imprint/Privacy