SDA (SDA flex)  7.2
Simulation of Diffusional Association
Functions/Subroutines
read_input_file.f90 File Reference

Functions/Subroutines

subroutine read_input_file (tab_protein, rcriteria, resid_time, o_complexes, trajectories, geom, p_restart, type_calc, param_timestep, param_probe, param_force_energy, param_analytic, param_metaldesolv, param_pmf, option_omp, timemax, rboost, nrun, nmax_overlap, account_occurency, old_rotation)
 Function to read the input file and allocate all proteins, grids, flexible and all other variables. More...
 
subroutine read_list_conf (plist_flex, input_param, counter_grid, filename, param_probe, param_analytic, tmp_exclud, center, rcriteria, type_calc, been_read, fixed_cent)
 Read a list of conformations for flexible solutes. More...
 
subroutine initialize_bits (o_record, tab_prot)
 Initialize bits. More...
 
subroutine initialize_movable_atoms (type_calc, tab_prot)
 Initialize arrays for movable atoms. More...
 
subroutine initialize_random_generator (type_calc, dseed)
 Initialize the random number generators. More...
 

Detailed Description

Version
{version 7.2.3 (2019)}

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


Main function to read sda input file

Function/Subroutine Documentation

◆ initialize_bits()

subroutine initialize_bits ( type ( record o_record,
type ( array_protein_type ), target  tab_prot 
)

Initialize bits.

Could be in mod_record array_prot is included now.

Todo:
to move to mod_record ?

Could use param_force_energy, but bits info needed for creating array for trajectory to check if possible to setup param_force_energy before any record is initialized

Initialize bits, depending on the loaded grids ( bit_energy ) and type type calculation ( bit_integer )

Parameters
o_record: instance of mod_record::record
tab_prot: instance of array_protein_type
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialize_movable_atoms()

subroutine initialize_movable_atoms ( type ( type_calculation type_calc,
type ( array_protein_type tab_prot 
)

Initialize arrays for movable atoms.

Initialize the pointers of array of movable atoms for proteins
Check if they can be common and allocated into mod_setofgrid::sogrid, or if a local copy is needed into each protein
By default, points to sogrid. If needed, make a local copy for each protein (sdamm, sda_flex protein 2). If possible delete the copy in sogrid ( sdamm )
Do not account for OpenMP here, will be done in copy_array_protein

Parameters
type_calc: instance of mod_type_calc::type_calculation "type_calculation"
tab_prot: instance of array_protein_type
Here is the caller graph for this function:

◆ initialize_random_generator()

subroutine initialize_random_generator ( type ( type_calculation type_calc,
real ( kind=8 )  dseed 
)

Initialize the random number generators.

If dseed is equal to zero, the time is used as a seed, otherwise use the provided dseed

Initialze both generator ( see RNG : Random Number Generator ):

  • without OpenMP : only one seed is used
  • with OpenMP : each thread N uses the seed dseed+N
Parameters
type_calc: instance of type_calculation
dseed: input seed for RNG, if 0 use the time
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_input_file()

subroutine read_input_file ( type ( array_protein_type tab_protein,
type ( react_criter ), target  rcriteria,
type ( residence_time resid_time,
type ( record o_complexes,
type ( record trajectories,
type ( geometry geom,
type ( record ), pointer  p_restart,
type ( type_calculation type_calc,
type ( parameter_timestep param_timestep,
type ( probe_type param_probe,
type ( type_force_energy param_force_energy,
type ( parameter_analytic param_analytic,
type ( parameter_metaldesolv param_metaldesolv,
type ( parameter_pmf param_pmf,
type ( type_option_omp option_omp,
real ( kind=8 ), intent(out)  timemax,
real ( kind = 4 )  rboost,
integer  nrun,
integer  nmax_overlap,
logical  account_occurency,
logical, intent(out)  old_rotation 
)

Function to read the input file and allocate all proteins, grids, flexible and all other variables.

Print summary informations about parameters

Simple variable read by read_single_parameter :

  • number of runs, timemax, probe, scale factor for grids...

Use of groups/block for complexe modules, read by read_group :

  • No maximum predefined number of Solute_Grid groups
  • No maximum predefined Flexible list and proteins
  • Could be extended to more than one reaction criteria ??

Inconvenient:

  • Order of some GROUPs are important
Parameters
tab_protein: instance of array_protein_type
rcriteria: instance of react_criter
resid_time: instance of residence_time
o_complexes,trajectories,p_restart: instance of record
geom: instance of geometry
type_calc: instance of type_calculation
param_timestep: instance of parameter_timestep
param_probe: instance of probe_type
param_force_energy: instance of type_force_energy
param_analytic: definition of analytical interaction, instance of parameter_analytic
param_metaldesolv: instance of parameter_metaldesolv
option_omp: used for OpenMP
timemax: maximum time for each trajectory, 0 is infinite
rboost: distance when boost
nrun: number of trajectories
nmax_overlap: number of overlap before applying a boost
account_occurency: if occurences must be accounted
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_list_conf()

subroutine read_list_conf ( type ( list_flex_type plist_flex,
type ( param_list input_param,
integer  counter_grid,
character ( 128 )  filename,
type ( probe_type param_probe,
type ( parameter_analytic param_analytic,
logical  tmp_exclud,
real ( kind=8 ), dimension(3)  center,
type ( react_criter rcriteria,
type ( type_calculation type_calc,
logical  been_read,
real ( kind=8 ), dimension( 3 ), intent(in)  fixed_cent 
)

Read a list of conformations for flexible solutes.

Read an external input file
Two methods implemented :

  • EXPLICIT : give the names of all grids explicitely
  • IMPLICIT : use an automatic numerotation ( eg : pbd1.pdb, pdb2.pdb,... ), need additional parameters

If not flexible return the correct name of the grid ( used if fparallel )

Parameters
plist_flex: instance of list_flex_type
input_param: instance of param_list
counter_grid: use as identifier for grids
filename: name of external file with the list of conformations
param_probe: instance of probe_type
param_analytic: structure parameter_analytic
tmp_exclud: if exclusion must be computed
center: center of solutes ?
rcriteria: instance of react_criter
type_calc: instance of type_calculation
been_read: used for checking misspelling
Here is the call graph for this function:
Here is the caller graph for this function:
Imprint/Privacy