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

Module to define one protein (position, orientation, is surface?)
. More...

Data Types

type  protein
 Define properties of one solute / protein. More...
 

Functions/Subroutines

subroutine new_protein (this, irot, isurf, iflex, iimg_chg, id, vert_excl, opt_imdesolv)
 Allocate one protein. More...
 
subroutine set_flexible (this, plist_flexible)
 Allocate a flexible and assign the plist_flexible. More...
 
subroutine allocate_array_mdesolv (this, param_metaldesolv)
 
subroutine delete_protein (this)
 Destructor, normally not called directly but by mod_array_protein add new variable local_copy, simpler and allow easier for flex. More...
 
subroutine update_rot_orientation (this, mrot)
 rotate orientation axis version with explicit terms, faster More...
 
subroutine update_rotation_protein (this, matrix_rot)
 Rotate the position of atoms/ charges / etc. version with explicit rotation, faster now use the pointer, either points to a sogrid copy, or rotate the local copy. More...
 
subroutine info_protein (this)
 Print info. More...
 
subroutine init_movable_array (this, sgrid, make_local_copy)
 Make the copy of the atomic position to protein if make_local_copy is true Otherwise points to sogrid why sgrid ? already associated ? Not if flex to check. More...
 
subroutine copy_movable_array (this, opt_mat)
 Copy movable atoms from sogrid arrays to a local copy in the protein. Needed in case of flexibility (keep a copy of the initial position), certianly sda_pmf (with option opt_mat)\ Check if local copy is true, otherwise do nothing
Optional opt_mat will load the initial position locally and apply the matrix of rotation.
. More...
 
subroutine select_charge (grid_number, prot1, prot2, nq1, nq2, p_charge1, p_charge2, p_charge_position1, p_charge_position2, dummy_array)
 First version, fastest version, do not check for asymetry.
. More...
 
subroutine select_charge2 (grid_number, prot1, prot2, nq1, nq2, p_charge1, p_charge2, p_charge_position1, p_charge_position2, dummy_array, opt_pgrid1, opt_pgrid2, opt_subgrid)
 Select correct pointer for charges ( and optionaly subgrids : only lennard-jones now ).
. More...
 
subroutine initialize_conformation (this, opt_conf_nb)
 Initialize the conformation in case of flexible.
replace change_conformation::init_conformation
. More...
 
subroutine update_visited_conf (this, timestep)
 update array statistics for sda_flex, use timestep rather than integer, may change More...
 
subroutine merge_visited_conf (this, to_merge)
 
subroutine print_visited_conf (this)
 

Detailed Description

Module to define one protein (position, orientation, is surface?)
.

Function/Subroutine Documentation

◆ allocate_array_mdesolv()

subroutine mod_protein::allocate_array_mdesolv ( type ( protein this,
type ( parameter_metaldesolv param_metaldesolv 
)
Here is the caller graph for this function:

◆ copy_movable_array()

subroutine mod_protein::copy_movable_array ( type ( protein this,
real(kind=8), dimension(3,3), optional  opt_mat 
)

Copy movable atoms from sogrid arrays to a local copy in the protein. Needed in case of flexibility (keep a copy of the initial position), certianly sda_pmf (with option opt_mat)\ Check if local copy is true, otherwise do nothing
Optional opt_mat will load the initial position locally and apply the matrix of rotation.
.

Precondition
To be called after init_movable_array
Here is the call graph for this function:
Here is the caller graph for this function:

◆ delete_protein()

subroutine mod_protein::delete_protein ( type ( protein this)

Destructor, normally not called directly but by mod_array_protein add new variable local_copy, simpler and allow easier for flex.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ info_protein()

subroutine mod_protein::info_protein ( type ( protein ), intent(in)  this)

Print info.

◆ init_movable_array()

subroutine mod_protein::init_movable_array ( type ( protein this,
type ( sogrid ), pointer  sgrid,
logical  make_local_copy 
)

Make the copy of the atomic position to protein if make_local_copy is true Otherwise points to sogrid why sgrid ? already associated ? Not if flex to check.

Here is the caller graph for this function:

◆ initialize_conformation()

subroutine mod_protein::initialize_conformation ( type ( protein this,
integer, optional  opt_conf_nb 
)

Initialize the conformation in case of flexible.
replace change_conformation::init_conformation
.

Parameters
this: instance of protein
opt_conf_nb:: optional, indicates which conformation. If zero chosses a random one
Here is the call graph for this function:

◆ merge_visited_conf()

subroutine mod_protein::merge_visited_conf ( type ( protein this,
type ( protein to_merge 
)

◆ new_protein()

subroutine mod_protein::new_protein ( type ( protein this,
logical, intent(in)  irot,
logical, intent(in)  isurf,
logical, intent(in)  iflex,
logical, intent(in)  iimg_chg,
integer, intent(in)  id,
logical, intent(in)  vert_excl,
logical, optional  opt_imdesolv 
)

Allocate one protein.

Parameters
this: instance of protein
irot: if the solute rotates
isurf: if the solute is a surface ( no rotation, no translation, different PBC )
iflex: if the solute is flexible
iimg_chg: if solute should compute image_charge forces and energies
id: identifier, depreceated
vert_excl: should decision on whether to make an excluded volume check depend on vertical height
opt_imdesolv: optional, global object mdesolv
Here is the caller graph for this function:

◆ print_visited_conf()

subroutine mod_protein::print_visited_conf ( type ( protein this)

◆ select_charge()

subroutine mod_protein::select_charge ( integer, intent(in)  grid_number,
type ( protein ), intent(in), target  prot1,
type ( protein ), intent(in), target  prot2,
integer, intent(out)  nq1,
integer, intent(out)  nq2,
real ( kind = 8 ), dimension ( : ), intent(out), pointer  p_charge1,
real ( kind = 8 ), dimension ( : ), intent(out), pointer  p_charge2,
real ( kind=8 ), dimension(:,:), intent(out), pointer  p_charge_position1,
real ( kind=8 ), dimension(:,:), intent(out), pointer  p_charge_position2,
real ( kind = 8 ), dimension ( : ), target  dummy_array 
)

First version, fastest version, do not check for asymetry.
.

Can be called by sdamm, because an additional check is done in compute_force/energy

Assign pointer to the correct grid / charge.

Parameters
grid_number: grid to select
prot1,prot2: protein
nq1,nq2: number of "charges"
p_charge1,p_charge2: pointer to the "charge" value
p_charge_position1,p_charge_position2: pointer to the "charge" positions
dummy_array: sned the dummy array adress, used by hydrophobic and lennard-jones
Here is the caller graph for this function:

◆ select_charge2()

subroutine mod_protein::select_charge2 ( integer, intent(in)  grid_number,
type ( protein ), intent(in), target  prot1,
type ( protein ), intent(in), target  prot2,
integer, intent(out)  nq1,
integer, intent(out)  nq2,
real ( kind = 8 ), dimension ( : ), intent(out), pointer  p_charge1,
real ( kind = 8 ), dimension ( : ), intent(out), pointer  p_charge2,
real ( kind=8 ), dimension(:,:), intent(out), pointer  p_charge_position1,
real ( kind=8 ), dimension(:,:), intent(out), pointer  p_charge_position2,
real ( kind = 8 ), dimension ( : ), target  dummy_array,
type ( type_grid ), intent(out), optional, pointer  opt_pgrid1,
type ( type_grid ), intent(out), optional, pointer  opt_pgrid2,
integer, intent(in), optional  opt_subgrid 
)

Select correct pointer for charges ( and optionaly subgrids : only lennard-jones now ).
.

Include lennard-jones grids
Check if associated, allow asymetric computation
Could replace select_charge, need benchmark.

Parameters
grid_number: grid to select
prot1,prot2: protein
nq1,nq2: number of "charges"
p_charge1,p_charge2: pointer to the "charge" value
p_charge_position1,p_charge_position2: pointer to the "charge" positions
dummy_array: sned the dummy array adress, used by hydrophobic and lennard-jones
opt_pgrid1,opt_pgrid2: optional, pointer to the grid, can be lennard-jones
opt_subgrid: optional, the subgrid number ( only lj now )
Here is the caller graph for this function:

◆ set_flexible()

subroutine mod_protein::set_flexible ( type ( protein this,
type ( list_flex_type plist_flexible 
)

Allocate a flexible and assign the plist_flexible.

Parameters
this: instance of protein
plist_flexible: pointer to a list_flex_type
Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_rot_orientation()

subroutine mod_protein::update_rot_orientation ( type ( protein this,
real ( kind=8 ), dimension(3,3), intent(in)  mrot 
)

rotate orientation axis version with explicit terms, faster

Here is the caller graph for this function:

◆ update_rotation_protein()

subroutine mod_protein::update_rotation_protein ( type ( protein this,
real ( kind=8 ), dimension(3,3), intent(in)  matrix_rot 
)

Rotate the position of atoms/ charges / etc. version with explicit rotation, faster now use the pointer, either points to a sogrid copy, or rotate the local copy.

Here is the caller graph for this function:

◆ update_visited_conf()

subroutine mod_protein::update_visited_conf ( type ( protein this,
real ( kind=4 )  timestep 
)

update array statistics for sda_flex, use timestep rather than integer, may change

Imprint/Privacy