Quasi-adiabatic custom potential¶
-
namespace mdk
-
struct QAContact
- #include <QuasiAdiabatic.hpp>
A struct detailing the present QA contact.
Public Types
-
enum Status
Enum for the status of a QA contact.
Values:
-
enumerator FORMING
The contact is forming.
-
enumerator BREAKING
The contact is breaking.
-
enumerator REMOVED
The contact has been removed (but since we don’t want to actually remove items from the vector, we just mark it this way.)
-
enumerator FORMING
-
enum Status
-
struct QAFreePair
- #include <QuasiAdiabatic.hpp>
A struct detailing a free pair that can form a QA contact.
Public Types
-
enum Status
Enum for the status of the free pair.
Values:
-
enumerator FREE
The pair is indeed free.
-
enumerator TAKEN
The pair has formed a QA contact, but since we don’t want to outright remove an element from a list, we only mark it as removed.
-
enumerator FREE
Public Members
-
int i1
Index of the first residue.
-
int i2
Index of the second residue.
-
Status status
Status of the pair.
-
enum Status
-
struct QADiff
- #include <QuasiAdiabatic.hpp>
A “diff” with respect to the QA contact. Since the formation of the contact requires there to be enough “stat slots”, and it’s possible that in an asynchronous execution more contacts would be formed while they shouldn’t, i.e. oversaturate the pool of stat slots, we separate the formation into a “virtual formation” which is done asynchronously and which adds a
QADiff
to a list of diffs, and a synchronous part which adds these diffs one after another, thus preventing an oversaturation.Public Functions
-
class QuasiAdiabatic : public mdk::NonlocalForce
- #include <QuasiAdiabatic.hpp>
The quasi-adiabatic potential.
Public Functions
-
virtual void bind(Simulation &simulation) override
Bind the force to the simulation. In particular, we fetch the
types
andchains
, load sidechain distances etc.- Parameters
simulation – Simulation to bind to.
-
virtual void asyncPart(Dynamics &dynamics) override
Asynchronous part of the computation. In particular, we compute forces between currently-existing pairs and check which of the free pairs could form contacts.
- Parameters
dynamics – Dynamics object to add potential energy and forces to.
-
virtual void syncPart(Dynamics &dynamics) override
Synchronous part of the computation. We sort the list of potentially-added contacts and try to add them sequentially.
- Parameters
dynamics – Dynamics object to add potential energy and forces to; here it should be unused unless we want to account for the forces and potential energy of the added contacts.
-
virtual void vlUpdateHook() override
Action to be performed when the Verlet list is updated. Here we update the lists of pairs in contact and free pairs, preserving the pairs that were in contact in the old list and is present in the new Verlet list.
Private Functions
-
virtual vl::Spec spec() const override
Generate a spec for the Verlet list.
- Returns
Generated spec for the Verlet list.
-
void computeNH()
Compute a list of vectors $n_i$ and $h_i$.
-
bool geometryPhase(vl::PairInfo const &p, QADiff &diff) const
Perform a geometry check between two residues during the formation pass.
- Parameters
p – Data regarding the pair currently processed.
diff – A
QADiff
structure to modify if the check is successful.
- Returns
true
if the geometry check has succeeded,false
otherwise.
Private Members
-
Types const *types = nullptr
Types of residues, used when selecting the sidechain-sidechain potential.
-
Chains const *chains = nullptr
Chains data structure, here used mostly to filter the residues that are at the ends of the residues (there $n_i$ and $h_i$ cannot be computed).
-
LennardJones bb_lj = LennardJones(5.0 * angstrom, 1.0 * eps)
Lennard-Jones potential for the backbone-backbone contacts.
-
LennardJones bs_lj = LennardJones(6.8 * angstrom, 1.0 * eps)
Lennard-Jones potential for the backbone-sidechain contacts.
-
SidechainLJ ss_ljs[AminoAcid::N][AminoAcid::N]
An array of sidechain L-J potentials, one for every pair of types of amino acids.
-
Vectors n
An array of vectors $n_i$, as defined in CPC14.pdf.
-
Vectors h
An array of vectors $h_i$, as defined in CPC14.pdf.
-
double formationMaxDistSq
A maximum distance (across all types of contacts) between residues for the formation to occur. Has no “physical” sense, is used chiefly for optimization.
-
double hr_abs_min = 0.92
Minimum value of <h_i, r_j> for the formation of a backbone contact (with backbone part for i).
-
double hh_abs_min = 0.75
Minimum value of <h_i, h_j> for the formation of a backbone-backbone contact to occur.
-
double nr_max = 0.5
Maximum value of <n_i, r_j> for the formation of a sidechain contact (with sidechain part for i).
-
double formationTolerance = 1.0
A scalar factor by which the “usual” minimum formation distance is multiplied.
-
double formationTime = 10.0 * tau
A timespan during which a formed QA contact “thermalizes”, i.e. reaches the full depth of the L-J potential. The increase of the depth is linear.
-
double breakingTolerance = 1.0
A scalar factor by which the “usual” maximum distance before breaking is multiplied.
-
double breakingTime = 10.0 * tau
Timespan during which a broken QA contact “dissipates”, i.e. loses its force. The decrease in the depth of the L-J potential is linear.
-
std::vector<QAContact> oldPairs
The list of old pairs, with which the new pairs are swapped. This is done in order to not have to allocate new memory each time a Verlet list is regenerated.
-
std::vector<QAFreePair> freePairs
A list of free pairs, i.e. ones which are not in a QA contact but which potentially could form one.
-
Stats *stats
A reference to
Stats
variable. Note: this is a non-const reference, and in particular it gets potentially modified in thesyncPart
.
-
std::mutex qaDiffsMutex
A mutex for the access to the
qaDiffs
list. We have decided not to employ a separate list for every thread since (we estimate) there shouldn’t be many conflicts to the access; nevertheless we wish to prevent any concurrent access were one to occur.
-
virtual void bind(Simulation &simulation) override
-
struct QAContact