Data Types | |
interface | change_conformation |
called if flexible More... | |
Functions/Subroutines | |
subroutine | mainloop_sdamm (tab_protein, rcriteria, o_complexes, trajectories, geom, resid_time, type_calc, param_timestep, param_force_energy, param_analytic, param_metaldesolv, timemax, old_rotation) |
Main function for sdamm. More... | |
subroutine | compute_rdf (resid_time, tab_protein, geom) |
Call during the simulation to compute the radial distribution function. 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
Loop in the case of multiple proteins
subroutine mainloop_sdamm::compute_rdf | ( | type ( residence_time ) | resid_time, |
type ( array_protein_type ), intent(in) | tab_protein, | ||
type ( geometry ) | geom | ||
) |
Call during the simulation to compute the radial distribution function.
Certainly not optimal
resid_time | : instance residence_time |
tab_protein | : instance of array_protein_type |
geom | : instance of geometry |
subroutine mainloop_sdamm | ( | tab_protein, | |
rcriteria, | |||
o_complexes, | |||
trajectories, | |||
geom, | |||
resid_time, | |||
type_calc, | |||
param_timestep, | |||
param_force_energy, | |||
param_analytic, | |||
param_metaldesolv, | |||
timemax, | |||
old_rotation | |||
) |
Main function for sdamm.
In SDAMM ( SDA with Multiple Molecules ) we simulate crowding effects
Many proteins are inside a box ( with or without PBC (periodic boundary conditions ) and are simulated for a definite time
The main output is a trajectory file
The case where a solute is a surface is treated slightly differently ( use 2 cutoffs, modify PBC )
Some, or all solutes may also be flexible
Only parts of the computation is parrallelised :
tab_protein | : instance of array_protein_type |
rcriteria | : reaction criteria |
o_complexes | : instance of record for complex output |
trajectories | : instance of record for trajectory output |
geom | : instance of geometry |
resid_time | : instance of residence_time, here used for radial distribution function |
type_calc | : instance of type_calculation |
param_timestep | : instance of parameter_timestep |
param_force_energy | : structure type_force_energy |
param_analytic | : instance of parameter_analytic |
param_metaldesolv | : instance of parameter_metaldesolv |
timemax | : length of the trajectory |