Calculating entropy of DNA strand (Solute)

Loading Data

  1. Load your data into a MDAnalysis Universe.

[1]:
import MDAnalysis as mda
# set the working dir to the root of repo inorder to use these path
tprfile = "data/md_A4_dna.tpr"
trrfile = "data/md_A4_dna_xf.trr"
u = mda.Universe(tprfile, trrfile)
  1. Trim the data to reduce analysis time

[2]:
from CodeEntropy.IO import MDAUniverseHelper as MDAHelper
# since this trajectory only contains the DNA strand
selection_string_pre_process = 'all'
start = 3
end = 40
step = 1

reduced_frame = MDAHelper.new_U_select_frame(u,  start, end, step)

reduced_atom = MDAHelper.new_U_select_atom(reduced_frame, selection_string_pre_process)


  1. parse the data into a CodeEntropy data object

[3]:
from CodeEntropy.ClassCollection import DataContainer as DC
dataContainer = DC.DataContainer(reduced_atom)
Number of atoms      : 254
Number of frames     : 37

Performing calculation

The total entropy for a system is taken as the sum of seven terms: \(S_{total} = S^{transvib}_{WM} + S^{rovib}_{WM} + S^{transvib}_{R} + S^{rovib}_{R} + S^{transvib}_{UA} + S^{rovib}_{UA} + S^{topo}_{UA}\)

Set Parameters

[4]:
from CodeEntropy.FunctionCollection import EntropyFunctions as EF
#Only the part of interest remained
selection_string = "all"
axis_list = ["C5'", "C4'", "C3'"]
outfile = None
moutFile = None
nmdFile = None
csv_out = None
tScale = 1.0
fScale = 1.0
temper = 300.0 #K
verbose = 3
thread = 4

Whole-molecule Level

\(S^{transvib}_{WM} + S^{rovib}_{WM}\)

[5]:
wm_entropyFF, wm_entropyTT = EF.compute_entropy_whole_molecule_level(
    arg_hostDataContainer = dataContainer,
    arg_outFile = outfile,
    arg_selector = selection_string,
    arg_moutFile = moutFile,
    arg_nmdFile = nmdFile,
    arg_fScale = fScale,
    arg_tScale = tScale,
    arg_temper = temper,
    arg_verbose = verbose
)
print(f"wm_entropyFF = {wm_entropyFF}")
print(f"wm_entropyTT = {wm_entropyTT}")
------------------------------------------------------------
          Hierarchy level. --> Whole molecule <--
------------------------------------------------------------
Total number of beads at the whole molecule level = 1
Assigning Translation and Rotation Axes @ whole molecule level-> Done
Updating Local coordinates-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (Whole mol level)             : 48.3941 J/mol/K
TT Entropy (Whole mol level)             : 43.5001 J/mol/K
wm_entropyFF = 48.39410817515473
wm_entropyTT = 43.50008156137956

Residue level

\(S^{transvib}_{R} + S^{rovib}_{R}\)

[6]:
res_entropyFF, res_entropyTT = EF.compute_entropy_residue_level(
    arg_hostDataContainer = dataContainer,
    arg_outFile = outfile,
    arg_selector = selection_string,
    arg_moutFile = moutFile,
    arg_nmdFile = nmdFile,
    arg_fScale = fScale,
    arg_tScale = tScale,
    arg_temper = temper,
    arg_verbose = verbose,
    arg_axis_list = axis_list,
)

print(f"res_entropyFF = {res_entropyFF}")
print(f"res_entropyTT = {res_entropyTT}")
------------------------------------------------------------
             Hierarchy level. --> Residues <--
------------------------------------------------------------
DA1
DA2
DA3
DA4
DT5
DT6
DT7
DT8
Total number of beads at the residue level = 8
Assigning Translation Axes @ residue level-> Done
Assigning Rotational Axes @ residue level->
    1     2     3     4     5
    6     7     8
Done
Updating Local forces-> Done
Updating Local coordinates-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (Residue level)               : 201.0312 J/mol/K
TT Entropy (Residue level)               : 311.2150 J/mol/K
res_entropyFF = 201.03116269045088
res_entropyTT = 311.2149733303239

United Atom Level

\(S^{transvib}_{UA} + S^{rovib}_{UA}\)

[7]:
UA_entropyFF, UA_entropyTT, res_df = EF.compute_entropy_UA_level(
    arg_hostDataContainer = dataContainer,
    arg_outFile = outfile,
    arg_selector = selection_string,
    arg_moutFile = moutFile,
    arg_nmdFile = nmdFile,
    arg_fScale = fScale,
    arg_tScale = tScale,
    arg_temper = temper,
    arg_verbose = verbose,
    arg_csv_out= csv_out,
    arg_axis_list = axis_list,
)
print(f"UA_entropyFF = {UA_entropyFF}")
print(f"UA_entropyTT = {UA_entropyTT}")
print("Per residue entropy:")
print(res_df)
------------------------------------------------------------
            Hierarchy level. --> United Atom <--
------------------------------------------------------------
Working on resid : DA1
Total number of UA beads in residue DA1 : 18
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (54,54) -> (21, 21)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DA1)                  : 44.5105 J/mol/K
TT Entropy (UA for DA1)                  : 22.6489 J/mol/K



Working on resid : DA2
Total number of UA beads in residue DA2 : 21
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (63,63) -> (19, 19)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DA2)                  : 34.7319 J/mol/K
TT Entropy (UA for DA2)                  : 20.6152 J/mol/K



Working on resid : DA3
Total number of UA beads in residue DA3 : 21
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (63,63) -> (19, 19)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DA3)                  : 36.1097 J/mol/K
TT Entropy (UA for DA3)                  : 23.6253 J/mol/K



Working on resid : DA4
Total number of UA beads in residue DA4 : 21
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (63,63) -> (21, 21)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DA4)                  : 38.2059 J/mol/K
TT Entropy (UA for DA4)                  : 23.6248 J/mol/K



Working on resid : DT5
Total number of UA beads in residue DT5 : 17
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (51,51) -> (21, 21)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DT5)                  : 57.6747 J/mol/K
TT Entropy (UA for DT5)                  : 36.7715 J/mol/K



Working on resid : DT6
Total number of UA beads in residue DT6 : 20
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (60,60) -> (19, 19)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DT6)                  : 47.1265 J/mol/K
TT Entropy (UA for DT6)                  : 26.3445 J/mol/K



Working on resid : DT7
Total number of UA beads in residue DT7 : 20
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (60,60) -> (19, 19)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DT7)                  : 44.4613 J/mol/K
TT Entropy (UA for DT7)                  : 27.3633 J/mol/K



Working on resid : DT8
Total number of UA beads in residue DT8 : 20
Assigning Translation Axes at the UA level-> Done
Assigning Rotational Axes at the UA level-> Done
Updating Local forces-> Done
Updating Local torques-> Done
Weighting forces and torques-> Done
Updating the submatrices ...
Finished updating submatrix FF
Finished updating submatrix TT
Done
Generating Quadrants-> Done
A shape change has occured (60,60) -> (21, 21)
Diagonalizing-> Done
Changing the units of eigen values to SI units-> Done
Sorting spectrum in ascending order of frequencies-> Done
Entropy values:
FF Entropy (UA for DT8)                  : 39.4173 J/mol/K
TT Entropy (UA for DT8)                  : 33.6917 J/mol/K



------------------------------------------------------------
Total Entropy FF (UA level) :         342.238 J/mol/K
Total Entropy TT (UA level) :         214.685 J/mol/K
------------------------------------------------------------
UA_entropyFF = 342.2378738899839
UA_entropyTT = 214.68518956812613
Per residue entropy:
  RESNAME  RESID  FF_ENTROPY(J/mol/K)  TT_ENTROPY(J/mol/K)
0      DA      1            44.510450            22.648880
1      DA      2            34.731911            20.615160
2      DA      3            36.109735            23.625295
3      DA      4            38.205923            23.624781
4      DT      5            57.674748            36.771516
5      DT      6            47.126494            26.344506
6      DT      7            44.461274            27.363316
7      DT      8            39.417339            33.691735

United Atom Level Multi process

Use Multiple thread to speed up operation

[8]:
UA_entropyFF, UA_entropyTT, res_df = EF.compute_entropy_UA_level_multiprocess(
    arg_hostDataContainer = dataContainer,
    arg_outFile = outfile,
    arg_selector = selection_string,
    arg_moutFile = moutFile,
    arg_nmdFile = nmdFile,
    arg_fScale = fScale,
    arg_tScale = tScale,
    arg_temper = temper,
    arg_verbose = verbose,
    arg_csv_out= csv_out,
    arg_axis_list = axis_list,
    arg_thread= thread,
)
print(f"UA_entropyFF = {UA_entropyFF}")
print(f"UA_entropyTT = {UA_entropyTT}")
print("Per residue entropy:")
print(res_df)
Finished updating submatrix FF
Finished updating submatrix TT
A shape change has occured (54,54) -> (21, 21)



Finished updating submatrix FF
Finished updating submatrix TT
A shape change has occured (63,63) -> (19, 19)



Finished updating submatrix FF
Finished updating submatrix FF
Finished updating submatrix TT
Finished updating submatrix TT
A shape change has occured (63,63) -> (21, 21)
A shape change has occured (63,63) -> (19, 19)






Finished updating submatrix FF
Finished updating submatrix TT
A shape change has occured (51,51) -> (21, 21)



Finished updating submatrix FF
Finished updating submatrix TT
A shape change has occured (60,60) -> (19, 19)



Finished updating submatrix FF
Finished updating submatrix FF
Finished updating submatrix TT
Finished updating submatrix TT
A shape change has occured (60,60) -> (21, 21)
A shape change has occured (60,60) -> (19, 19)






  RESNAME  RESID  FF_ENTROPY(J/mol/K)  TT_ENTROPY(J/mol/K)
0      DA      1            44.510450            22.648880
1      DA      2            34.731911            20.615160
2      DA      3            36.109735            23.625295
3      DA      4            38.205923            23.624781
4      DT      5            57.674748            36.771516
5      DT      6            47.126494            26.344506
6      DT      7            44.461274            27.363316
7      DT      8            39.417339            33.691735
UA_entropyFF = 342.23787388998386
UA_entropyTT = 214.6851895681262
Per residue entropy:
  RESNAME  RESID  FF_ENTROPY(J/mol/K)  TT_ENTROPY(J/mol/K)
0      DA      1            44.510450            22.648880
1      DA      2            34.731911            20.615160
2      DA      3            36.109735            23.625295
3      DA      4            38.205923            23.624781
4      DT      5            57.674748            36.771516
5      DT      6            47.126494            26.344506
6      DT      7            44.461274            27.363316
7      DT      8            39.417339            33.691735

Topographical entropy

\(S^{topo}_{UA}\)

[9]:
result_entropyAEM = EF.compute_topographical_entropy_AEM(
    arg_hostDataContainer = dataContainer,
    arg_selector = selection_string,
    arg_outFile = outfile,
    arg_verbose = verbose
)

print(f"result_entropyAEM = {result_entropyAEM}")
------------------------------------------------------------
Topographical entropy of residue side chains
computed using all the dihedrals with AEM method
------------------------------------------------------------
----------Working on resid : 1 (DA)----------
Found 18 exclusive dihedrals in residue DA1
Found 18 dihedrals which collectively acquire 8 unique conformers
Residue Topographical Entropy from AEM (DA 1) : 13.4948
----------Working on resid : 2 (DA)----------
Found 20 exclusive dihedrals in residue DA2
Found 20 dihedrals which collectively acquire 4 unique conformers
Residue Topographical Entropy from AEM (DA 2) : 7.6987
----------Working on resid : 3 (DA)----------
Found 20 exclusive dihedrals in residue DA3
Found 20 dihedrals which collectively acquire 2 unique conformers
Residue Topographical Entropy from AEM (DA 3) : 2.8480
----------Working on resid : 4 (DA)----------
Found 19 exclusive dihedrals in residue DA4
Found 19 dihedrals which collectively acquire 2 unique conformers
Residue Topographical Entropy from AEM (DA 4) : 1.0331
----------Working on resid : 5 (DT)----------
Found 14 exclusive dihedrals in residue DT5
Found 14 dihedrals which collectively acquire 4 unique conformers
Residue Topographical Entropy from AEM (DT 5) : 7.8987
----------Working on resid : 6 (DT)----------
Found 16 exclusive dihedrals in residue DT6
Found 16 dihedrals which collectively acquire 2 unique conformers
Residue Topographical Entropy from AEM (DT 6) : 1.0331
----------Working on resid : 7 (DT)----------
Found 16 exclusive dihedrals in residue DT7
Found 16 dihedrals which collectively acquire 2 unique conformers
Residue Topographical Entropy from AEM (DT 7) : 1.0331
----------Working on resid : 8 (DT)----------
Found 15 exclusive dihedrals in residue DT8
Found 15 dihedrals which collectively acquire 1 unique conformers
Residue Topographical Entropy from AEM (DT 8) : -0.0000
------------------------------------------------------------
Total Topog. Entropy (AEM)               :          35.040
------------------------------------------------------------
result_entropyAEM = 35.039528380573834

Total Entropy

[10]:
total = wm_entropyFF + wm_entropyTT + res_entropyFF + res_entropyTT + UA_entropyFF + UA_entropyTT + result_entropyAEM
print(f"Total Entropy = {total} J/mol/K")
Total Entropy = 1196.1029175959927 J/mol/K
[ ]: