Potential restraints

Potential restraints are a family of EnergyFunctionComponent instances that apply any given potential to a selection of Atom instances, optionally further masked or mapped. In other words, the application of a potential restraint to a Pose is characterized by a series of fully customizable steps to generate new and/or improved Potential restraints:

  1. (Optional) Select a subset of Atom instances to apply the potential;
  2. Define the potential function to be applied (See Available potentials);
  3. Apply the potential mask to the selected Atom instances;
  4. (Optional) Multiply a mask/map to the calculated energy and forces matrices;
  5. Sum the calculated energy and forces matrices.
Note:

Most of these EnergyFunctionComponent instances are specific for a certain type of molecules (such as Peptides). The following types and methods constitute the backbone for all Potential restraints, but more specific applications can be found in other modules of ProtoSyn.

ProtoSyn Potential Restraint

Figure 1 | A diagram representation of the application of a potential restraint, in ProtoSyn. A given potential function (See Available potentials) is applied to a distance matrix (See Distance matrix calculation) calculated from the distances between all selected Atom instances. The resulting 2D energy matrix includes the energy values for each pair of interacting Atom instances. An optional secondary output from the application of the potential function is a 3D forces matrix, with the forces felt on each atom from the interaction with its pair. Both these matrices can optionally suffer another step: the application of a mask or map. On one hand, a Mask, similarly to its application in AbstractSelections, turns on and off certain pairs of interacting atoms. As an example, a mask allows the user to ignore same-atom energy artifacts using a get_diagonal_mask or ignore interaction between atoms of the same residue using an get_intra_residue_mask. On the other hand, a map in a 2D matrix (with the same size as the 2D energy matrix) with weights and biases. When the two matrices are multiplied together, certain interaction strengths are increased while others are reduced. An example would be the application of a contact map, where each weight corresponds to the certainty factor for the distance between each of the considered Atom pairs. As a final step both the 2D energy matrix and the 3D forces matrix (optionally after the application of the mask/map) are summed (in both axis, in the case of energy calculation, and in one axis, in the case of forces calculation).

Note:

The above paradigm for application of potential restraints requires, for the correct calculation of forces, the full_distance_matrix to be considered. This causes the energy value calculated to be double the expected (not considering the application on any non diagonally symmetrical mask), as both the top and bottom triagonal matrices are considered. In most cases this is not an issue, as the energy values are compared between simulation frames in relation to eachother and therefore the scale of the value is not important. However, in specific applications, further adjustment of the potential and/or mask applied may be necessary.

The following section is subdivided in the following topics, for organizational purposes:

Base methods

The following methods are ubiquous to all Potential restraints:

ProtoSyn.Calculators.apply_potential!Function
apply_potential!([A::Type{<: ProtoSyn.AbstractAccelerationType}], pose::Pose, potential::Function, update_forces::Bool, [verlet_list::Union{VerletList, Nothing}], [selection::Union{Nothing, ProtoSyn.AbstractSelection}], [mask::MaskMap])

Apply the given potential to the provided Pose pose. If the update_forces flag is set to true, also calculate the forces acting on the system. The call to this function can be further customized by providing an optional VerletList, an AbstractSelection selection or a MaskMap map (this can be a Mask, a Matrix or a Function).

Examples

julia> ProtoSyn.Calculators.apply_potential!(ProtoSyn.SISD_0, pose, (d; kwargs...) -> 0.0, false, nothing, nothing, nothing)
(0.0, [0.0 0.0 … 0.0 0.0; 0.0 0.0 … 0.0 0.0; 0.0 0.0 … 0.0 0.0])
source

Available potentials

A list of available potentials in ProtoSyn is described bellow. See Creating custom potential functions for the required signature of new potential functions.

ProtoSyn.Calculators.get_flat_bottom_potentialFunction
get_flat_bottom_potential(;d1::T = 0.0, d2::T = 0.0, d3::T = Inf, d4::T = Inf) where {T <: AbstractFloat}

Return a flat-bottom potential function, using the specified distances. The potential is made up of 5 different sectors, each with the following functions:

\[f_{1}) \;\;\;\;\;\; e = m_{1} \cdot d + b_{1} \;\;\;\;\;\; \left \{ d < d_{1} \right \}\,\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[f_{2}) \;\;\;\;\;\; e = \left (d-d_{2} \right )^{2} \;\;\;\;\;\;\;\; \left \{ d_{1} \leqslant d < d_{2} \right \}\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[f_{3}) \;\;\;\;\;\; e = 0 \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\; \left \{ d_{2} \leqslant d \leqslant d_{3} \right \}\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[f_{4}) \;\;\;\;\;\; e = \left (d-d_{3} \right )^{2} \;\;\;\;\;\;\;\; \left \{ d_{4} < d \leqslant d_{4} \right \}\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[f_{5}) \;\;\;\;\;\; e = m_{2} \cdot d + b_{2} \;\;\;\;\;\; \left \{ d > d_{4} \right \}\,\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

Where

\[m_{1} = 2 \left ( d_{1}-d_{2} \right ) \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[b_{1} = f_{2}\left ( d_{1} \right ) - m_{1} \cdot d_{1} \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[m_{2} = 2\left ( d_{4} - d_{3} \right ) \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

\[b_{2} = f_{4}\left ( d_{4} \right ) - m_{2} \cdot d_{4} \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\]

The resulting function can be called with the following signature:

flat_bottom_potential(d::T; v::Opt{Vector{T}} = nothing) where {T <: AbstractFloat}

Return an energy value based on the provided distance d. If a vector v is also provided (optional), the flat-bottom restraint will also return the forces f1 and f2 (the forces felt on both ends of the vector v). The vector v should have length = 3, corresponding to the 3 dimensions of the distance between the two Atom instances (X, Y and Z). For more information on the flat-bottom potential, please read: https://www.ncbi.nlm.nih.gov/pmc/articles/PMC4692055/.

See also

apply_potential! calc_flat_bottom_restraint

Examples

julia> f = ProtoSyn.Calculators.get_flat_bottom_potential(d1 = 1.0, d2 = 2.0, d3 = 3.0, d4 = 4.0);

julia> f(2.5)
0.0

julia> f(1.73, v = (1.0, 1.0, 1.0))
(0.0729, (0.54, 0.54, 0.54), (-0.54, -0.54, -0.54))
source
ProtoSyn.Calculators.get_coulomb_potentialFunction
get_coulomb_potential(d::T; v::Opt{Tuple{T, T, T}} = nothing, kwargs...) where {T <: AbstractFloat}

Return a simple coulomb potential with charges included, as described in https://www.softschools.com/formulas/physics/potentialenergyelectrostaticpointparticles_formula/37/. The coulomb potential is measured between two charged point-like particles at distance d and charges kwargs[:qi] and kwargs[:qj], respectivelly. If a vector v between the two particles is provided, also calculates the resulting forces.

See also

get_bump_potential

Examples

julia> coulomb = ProtoSyn.Calculators.get_coulomb_potential()

julia> coulomb(2.0, qi = 1.0, qj = -1.0)
-0.5
source
ProtoSyn.Calculators.get_bump_potentialFunction
get_bump_potential(;[c::T = 1.0], [r::T = 1.0]) where {T <: AbstractFloat}

Return a bump potential without charges included, as described in https://math.stackexchange.com/a/3236066. The bump function is centered around c with a radius r.

See also

get_bump_potential_charges

Examples

julia> ProtoSyn.Calculators.get_bump_potential(c = 1.0, r = 1.0)
(::ProtoSyn.Calculators.var"#bump_potential#25"{ProtoSyn.Calculators.var"#bump_potential#24#26"{Float64, Float64, Float64}}) (generic function with 1 method)
source
ProtoSyn.Calculators.get_bump_potential_chargesFunction
get_bump_potential_charges(;[c::T = 1.0], [r::T = 1.0]) where {T <: AbstractFloat}

Return a bump potential with charges included, as described in https://math.stackexchange.com/a/3236066. The bump function is centered around c with a radius r.

See also

get_bump_potential

Examples

julia> ProtoSyn.Calculators.get_bump_potential_charges(c = 1.0, r = 1.0)
(::ProtoSyn.Calculators.var"#bump_potential_charges#21"{ProtoSyn.Calculators.var"#bump_potential_charges#20#22"{Float64, Float64}}) (generic function with 1 method)
source
ProtoSyn.Calculators.get_harmonic_potentialFunction
get_harmonic_potential(a::T, b::T, c::T) where {T <: AbstractFloat}

Return a simple harmonic potential with the following equation:

U = ((d+a)²-b)*c

The harmonic potential is measured between two point-like particles at distance d (doesn't take charges into consideration). If a vector v between the two particles is provided, also calculates the resulting forces.

See also

get_bump_potential

Examples

julia> harmonic = ProtoSyn.Calculators.get_harmonic_potential(1.0, 1.0, 1.0)

julia> harmonic(1.0)
3.0
source

ProtoSyn Flat Bottom Potential

Figure 1 | Visual depiction of the various potentials made available with ProtoSyn (version 1.1).

Creating custom potential functions

Besides the currently Available potentials, users can freely create custom and/or improved potential functions. In order to incorporate these potentials in EnergyFunctionComponent instances (and therefore in EnergyFunction instances), the following signature must be followed.

energy, force1, force2 = new_potential(d::T; v::Opt{Tuple{T, T, T}} = nothing, kwargs...) where {T <: AbstractFloat}

Any function must, therefore, receive a distance value d::T between two Atom instances and return the corresponding energy felt. Optionally, it may receive a vector v::Vector{T}, which is the difference vector between the two Atom positions, in which case the potential should also return the two forces felt on each atom from this energetic interaction. The kwargs may contain extra information from the apply_potential! function. Currently, only the Atom charges are provided, with entries qi and qj. This may change in future versions of ProtoSyn.

The newly defined method can then be used in the apply_potential! function: for each pair of Atom instances in the full_distance_matrix this potential will be evaluated. The apply_potential! (with the new potential) can then be used in a new :calc field of any new EnergyFunctionComponent (with the correct signature, see Creating a custom EnergyFunctionComponent). However, as a general rule of good practice, an extra step is usually included: as exemplified by the get_flat_bottom_potential method, the potential here is actually a functor, returned everytime the get_flat_bottom_potential function is called. A set of settings act as the input of this function (the settings of the encompassing EnergyFunctionComponent), parametrizing the returned potential. Using this approach, the user has complete access to the EnergyFunctionComponent.settings field, which will dictate the parametrization of the newly developed potential. This can be useful, for example, in efforts to correctly parametrize an energetic contribution based on a set of data, as the used potential settings can be changed and optimized each step of the simulation, as a new potential is generated everytime. However, this generation should not negatively impact performance in any meaningfull way.

Available masks

A list of available masks/maps in ProtoSyn is described bellow. See Creating custom masks for the required signature of new potential function masks/maps, and consider using ProtoSyn.Calculators.show_available_masks method to list all available masks.

ProtoSyn.Calculators.show_available_masksFunction
show_available_masks([io::IO], [m::Module])

Prints all available masks for potential restraints Module m (defaults to ProtoSyn.Calculators) to the given IO io (defaults to stdout). Recursivelly searches any inner Module.

Examples

julia> ProtoSyn.Calculators.show_available_masks()
+-------------------------------------------------------+
| Index | Mask function                                 |
+-------------------------------------------------------+
| 1     | get_bonded_mask                               |
| 2     | get_diagonal_mask                             |
+-------------------------------------------------------+
  └── Consider using the `?` menu to learn more about each mask function.
source
ProtoSyn.Calculators.get_intra_residue_maskFunction
get_intra_residue_mask(pose::Pose, [selection::Opt{AbstractSelection}])

For all the Atom instances in the provided AbstractSelection selection (N), return a 2D N x N Mask with all the Atom instances of the given Pose pose not in the same residue selected.

Note:

This function is rather heavy and has low performance. If no design effort is being made (where the sequence changes), the resulting Mask from this function can and should be re-used (only calculated once). If, for a specific application, the AbstractSelection selection remains constant but the Mask needs to be re-calculated (for example, because there was a design/mutation step, use the functor resulting from get_intra_residue_mask).

See also

show_available_masks

Examples

julia> ProtoSyn.Calculators.get_intra_residue_mask(pose, !an"^CA$|^N$|^C$|^H$|^O$"r)
ProtoSyn.Mask
 ├── Type: Atom
 ├── Size: (1140, 1140)
 ├── Count: 1279946 / 1299600
 └── Content: [0 0 … 1 1; 0 0 … 1 1; … ; 1 1 … 0 0; 1 1 … 0 0]
source
ProtoSyn.Calculators.get_bonded_maskFunction
get_bonded_mask(pose::Pose, [selection::Opt{AbstractSelection}])

For all the atoms in the provided AbstractSelection selection (N), return a 2D N x N Mask: for each Atom instance of the given Pose pose mask out all other bonded Atom instance.

See also

show_available_masks

Examples

julia> ProtoSyn.Calculators.get_bonded_mask(pose, an"CA")
ProtoSyn.Mask
 ├── Type: Atom
 ├── Size: (73, 73)
 ├── Count: 5256 / 5329
 └── Content: [0 1 … 1 1; 1 0 … 1 1; … ; 1 1 … 0 1; 1 1 … 1 0]
source
ProtoSyn.Calculators.get_diagonal_maskFunction
get_diagonal_mask(pose::Pose, [selection::Opt{AbstractSelection}])

For all the atoms in the provided AbstractSelection selection (N), return a 2D N x N Mask with all the Atom instances of the given Pose pose not in the natural diagonal selected (i.e. ignores same Atom interaction artifacts).

See also

show_available_masks

Examples

julia> ProtoSyn.Calculators.get_diagonal_mask(pose, an"CA")
ProtoSyn.Mask
 ├── Type: Atom
 ├── Size: (73, 73)
 ├── Count: 5256 / 5329
 └── Content: [0 1 … 1 1; 1 0 … 1 1; … ; 1 1 … 0 1; 1 1 … 1 0]
source
ProtoSyn.Calculators.get_upper_triangular_matrix_maskFunction
get_upper_triangular_matrix_mask(pose::Pose, [selection::Opt{AbstractSelection}])

For all the atoms in the provided AbstractSelection selection (N), return a 2D N x N Matrix{T} with the bottom triangular matrix set to 0.0 (including diagonal) and upper triangular matrix set to 1.0.

See also

show_available_masks

Examples

julia> ProtoSyn.Calculators.get_upper_triangular_matrix_mask(pose, an"CA")
73×73 Matrix{Float64}:
 (...)
source
ProtoSyn.Calculators.get_upper_triangular_matrix_inversed_maskFunction
get_upper_triangular_matrix_mask(pose::Pose, [selection::Opt{AbstractSelection}])

For all the atoms in the provided AbstractSelection selection (N), return a 2D N x N Matrix{T} with the bottom triangular matrix set to 0.0 (including diagonal) and upper triangular matrix set to -1.0.

See also

show_available_masks

Examples

julia> ProtoSyn.Calculators.get_upper_triangular_matrix_inversed_mask(pose, an"CA")
73×73 Matrix{Float64}:
 (...)
source
ProtoSyn.Calculators.load_mapFunction
load_map([::Type{T}], filename::String) where {T <: AbstractFloat}

Load the map in the filename file (i.e. Contact Map). The file should be in PFRMAT RR format (See: https://predictioncenter.org/casp13/index.cgi?page=format#RR). Returns an N x N map of the found weights, with pairs not identified in the file set to 0.0 (N is the maximum indentifier found on the file. As an example, it might be the case where a peptide has 74 residues, but no pair with residue 74 is found on the file, the maximum identifier found might be 72, for example. In this case, the resulting map will have size 72 x 72. In order to ensure the loaded map size matches the underlying peptide size, consider adding an entry of 0.0 on the map file, with the correct maximum identifier). Note: If no optional type T is provided, will use ProtoSyn.Units.defaultFloat.

Examples

julia> cmap = ProtoSyn.Calculators.load_map("contact_map_example.txt")
73×73 Array{Float64,2}:
 (...)
source

ProtoSyn Available Masks

Figure 2 | Some of the available Mask instances in ProtoSyn. A get_intra_residue_mask de-selects all atoms that belong to the same Residue instance (only accounting for inter-Residue interactions). A get_diagonal_mask simply de-selects the same Atom instance in a 2D full_distance_matrix. Finally, a contact map is a loaded map (using the load_map method) that asserts a 2D Matrix of weights or biases to be multiplied by the resulting energy and forces matrices in apply_potential!. As an example, this can be the likelihood of two Atom instances having a distance shorter can D Angstrom (Å).

Creating custom masks

In addition to the Available masks, users may wish to create and add custom Mask instances to employ in the apply_potential! method. This method can receive a 2D Mask instance directly, as well as a Matrix{T} or a Function instance. Each case is discussed bellow:

  • 2D Mask : A static mask, calculated only once. Energy and forces contributions will only be considered for the selected Atom instances (make sure the size of the Mask matches the number of selected Atom instances in the apply_potential! call). As an example, see get_intra_residue_mask. The syntax to create a new ProtoSyn Mask is as follows. For more details, see the Masks section.
julia> 2d_mask = Mask{Atom}(BitArray(falses(4, 4)))ERROR: syntax: "2" is not a valid function argument name around REPL[1]:1
  • $Matrix{T}$ : A static map, calculated only once. Energy and forces contibutions will be multiplied by the Matrix{T} (make sure that the Matrix{T} size matches the number of selected Atom instances in the apply_potential! call). As an example, see load_map. The syntax to create a new Matrix{T} is as follows (as an example of a random map).
julia> 2d_map = rand(4, 4)ERROR: syntax: "2" is not a valid function argument name around REPL[1]:1
  • Function : A dynamic mask or map, calculated every step/call of apply_potential!. This function should have the following simple signature, receiving a Pose pose as input and returning either a Mask or a Matrix{T}.
mask_or_map(pose::Pose)

Altough potentially sacrificing the performance of the code, this allows for new Mask or a Matrix{T} instances to be calculated every step of a simulation, reflecting new changes, such as mutations or design efforts, for example.

Note:

A common practice in building new Function instances that return Mask or a Matrix{T} instances for application in the apply_potential! method is to encompass this Function as a functor of another Function. As an example, see the get_intra_residue_mask method. In this cases, the encompassing Function receives the AbstractSelection selection, making sure that the returned Function selects Atom instances only within this selection.

When a dynamic mask or map is used, performance penalties are to be expected. If a user whises to fixate a dynamic mask into a static mask (on the context of an EnergyFunctionComponent or EnergyFunction), the following companion methods are made available:

ProtoSyn.Calculators.fixate_mask!
ProtoSyn.Calculators.fixate_masks!

Available Potential Restraint EnergyFunctionComponents

Making use of all the above, the following default EnergyFunctionComponent instances are defined and made available:

ProtoSyn.Calculators.Restraints.calc_flat_bottom_restraintFunction
calc_flat_bottom_restraint([::Type{A}], pose::Pose, update_forces::Bool; d1::T = 0.0, d2::T = 0.0, d3::T = Inf, d4::T = Inf, selection::Opt{AbstractSelection} = nothing, mask::MaskMap = nothing) where {A <: ProtoSyn.AbstractAccelerationType, T <: AbstractFloat}

Apply a flat bottom potential to a given Pose pose. The potential is iteratively obtained each call using the get_flat_bottom_potential method (See Available potentials), by providing the given d1::T (default: 0.0), d2::T (default = 0.0), d3::T (default = Inf) and d4::T (default = Inf) settings as the flat bottom potential distances. This potential then applied to the Pose pose (via the apply_potential! method), optionally on a subset of Atom instances given by the AbstractSelection selection and optionally multiplied by a mask. This mask can be a Mask, a Matrix{T} or a Function, in which case it should be a functor (return a Function) (For the correct signature of this Function mask, see Creating custom masks). These 3 options are named MaskMap for a simplicity of organization only. Return the total energy of the system and matrix of forces felt on each atom. Note that the calculation acceleration type can be set by providing an option parameter A Type{<: ProtoSyn.AbstractAccelerationType}. If not provided, the default ProtoSyn.acceleration.active will be used instead.

calc_flat_bottom_restraint!([::Type{A}], pose::Pose, update_forces::Bool; d1::T = 0.0, d2::T = 0.0, d3::T = Inf, d4::T = Inf, selection::Opt{AbstractSelection} = nothing, mask::MaskMap = nothing) where {A <: ProtoSyn.AbstractAccelerationType, T <: AbstractFloat}

Apply a flat bottom potential to a given Pose pose (see above). Also apply any energy and forces changes directly to the Pose pose.

Examples

julia> ProtoSyn.Calculators.Restraints.calc_flat_bottom_restraint(pose, true)
(0.0, [0.0 0.0 … 0.0 0.0; 0.0 0.0 … 0.0 0.0; 0.0 0.0 … 0.0 0.0])

julia> ProtoSyn.Calculators.Restraints.calc_flat_bottom_restraint(pose, false, d1 = 10.0, d2 = 12.0)
(556449.1936070402, [-711.7603616347209 -630.2662235401388 … 995.0284325254745 1153.572133762037; -419.1275359380875 -548.0506257124055 … 286.5285847489888 92.16862928705675; 6.007398880372552 8.2409631821887 … -99.38257889245355 -92.37110004070036])    ```
source
ProtoSyn.Calculators.Restraints.get_default_all_atom_clash_restraintFunction
get_default_all_atom_clash_restraint(;[α::T = 1.0]) where {T <: AbstractFloat}

Return the all-atom clash EnergyFunctionComponent. α sets the component weight (on an EnergyFunction instance). This component employs the calc_flat_bottom_restraint method, therefore defining a Pose energy based on a flat-bottom potential function applied to all atom-pairs in the system (N² complexity). By default, this EnergyFunctionComponent potential sets :d1 and :d2 of the flat-bottom potential to be 1.0 and 2.0, and masks out bonded atom-pairs.

See also

get_default_bond_distance_restraint

Settings

  • d1::T - The :d1 distance in the flat-bottom potential;
  • d2::T - The :d2 distance in the flat-bottom potential;
  • d3::T - The :d3 distance in the flat-bottom potential;
  • d4::T - The :d4 distance in the flat-bottom potential;
  • mask::MaskMap - The Mask, Matrix{T} or Function (see Creating custom masks) that masks out our multiplied by a set of pre-defined Atom instances;

Examples

julia> ProtoSyn.Calculators.Restraints.get_default_all_atom_clash_restraint()
🞧  Energy Function Component:
+---------------------------------------------------+
| Name           | All_Atom_Clash_Rest              |
| Alpha (α)      | 1.0                              |
| Update forces  | true                             |
| Calculator     | calc_flat_bottom_restraint       |
+---------------------------------------------------+
 |    +----------------------------------------------------------------------------------+
 ├──  ● Settings                      | Value                                            |
 |    +----------------------------------------------------------------------------------+
 |    | d4                            | Inf                                              |
 |    | d2                            | 2.0                                              |
 |    | mask                          | get_bonded_mask                                  |
 |    | d1                            | 1.0                                              |
 |    | d3                            | Inf                                              |
 |    +----------------------------------------------------------------------------------+
 |    
 └──  ○  Selection: nothing
source