SDA (SDA flex)  7.2
Simulation of Diffusional Association
Data Types | Functions/Subroutines
MainLoop_sda_2proteins_omp.f90 File Reference

Data Types

interface  change_conformation
 called if flexible More...
 

Functions/Subroutines

subroutine mainloop_sda_2proteins_omp (tab_protein_unique, rcriteria, o_complexes, trajectories, geom, p_restart, resid_time, type_calc, param_timestep, param_probe, param_force_energy, param_analytic, param_metaldesolv, option_omp, nmax_overlap, rboost, timemax, nrun, old_rotation)
 MainLoop_sda_2proteins_omp. More...
 
subroutine convert_to_record_format (array_energy, array_energy_1d, max_bit_energy)
 Make format compatible to call add_record.
. 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


Parralel version of sda 2proteins with OpenMP

Main loop for BD with 2 solutes
Used for normal sda or sda_koff
Can compute association rate, docking or electron transfert calculation with 2 proteins,
Residence time, first passage time

Function/Subroutine Documentation

◆ convert_to_record_format()

subroutine mainloop_sda_2proteins_omp::convert_to_record_format ( real(kind=8), dimension(:,:), intent(in)  array_energy,
real(kind=8), dimension(:,:), intent(out)  array_energy_1d,
integer, intent(in)  max_bit_energy 
)

Make format compatible to call add_record.
.

could be in record, and get bit information from there. need a parameter, maximum max_bit_energy to limit the loop of the array.

Here is the caller graph for this function:

◆ mainloop_sda_2proteins_omp()

subroutine mainloop_sda_2proteins_omp (   tab_protein_unique,
  rcriteria,
  o_complexes,
  trajectories,
  geom,
  p_restart,
  resid_time,
  type_calc,
  param_timestep,
  param_probe,
  param_force_energy,
  param_analytic,
  param_metaldesolv,
  option_omp,
  nmax_overlap,
  rboost,
  timemax,
  nrun,
  old_rotation 
)

MainLoop_sda_2proteins_omp.

Main function for 2 solutes computation

Parameters
tab_protein_unique: contains all proteins, instance of array_protein_type
rcriteria: instance of react_criter
o_complexes: for writing complexes
trajectories: for writing trajectory
geom: geometry instance
p_restart: pointer on a restart record
resid_time: instance of residence_time
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: structure parameter_analytic
param_metaldesolv: instance of parameter_metaldesolv
option_omp: option for openmp
nmax_overlap: could with other parameter
rboost: distance of boost
timemax: maximum time for a trajectory, 0 means infinite
nrun: number of runs
old_rotationstates whether the old rotation should be used or the one via Rodrigues formula
Here is the caller graph for this function:
Imprint/Privacy