Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fixing coordnum #25

Merged
merged 7 commits into from
Dec 5, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 7 additions & 2 deletions src/mddf.jl
Original file line number Diff line number Diff line change
Expand Up @@ -176,13 +176,13 @@ function mddf(trajectory::Trajectory, options::Options = Options(); coordination
# Read frame coordinates
lock(read_lock) do
# skip frames if stride > 1
while (iframe + 1) % options.stride != 0
while (iframe - options.firstframe + 1) % options.stride != 0
nextframe!(trajectory)
iframe += 1
end
# Read frame for computing
iframe += 1
nextframe!(trajectory)
iframe += 1
# The solute coordinates must be read in intermediate arrays, because the
# solute molecules will be considered one at a time in the computation of the
# minimum distances
Expand Down Expand Up @@ -494,6 +494,11 @@ end
R2 = mddf(traj2, Options(lastframe=2, frame_weights=[0.5, 0.25]))
@test R2.md_count ≈ R1.md_count
@test R2.volume.total ≈ R1.volume.total
# Varying weights with stride
R1 = mddf(traj1, Options(firstframe=2, lastframe=3))
R2 = mddf(traj1, Options(lastframe=3, stride=2, frame_weights=[1.0, 100.0, 1.0]))
@test R2.md_count ≈ R1.md_count
@test R2.volume.total ≈ R1.volume.total

end

Expand Down
6 changes: 6 additions & 0 deletions src/results.jl
Original file line number Diff line number Diff line change
Expand Up @@ -402,6 +402,7 @@ function finalresults!(R::Result, options::Options, trajectory::Trajectory)
#
# Computing the distribution functions and KB integrals, from the MDDF and from the RDF
#
warn = false
for ibin = 1:R.nbins
# For the MDDF
if R.md_count_random[ibin] > 0.0
Expand All @@ -412,6 +413,11 @@ function finalresults!(R::Result, options::Options, trajectory::Trajectory)
for j = 1:trajectory.solvent.natomspermol
R.solvent_atom[ibin, j] = R.solvent_atom[ibin, j] / R.md_count_random[ibin]
end
else
if !warn && !options.silent
@warn "Ideal-gas histogram bins with zero samples. Increase n_random_samples or trajectory length."
warn = true
end
end
if ibin == 1
R.coordination_number[ibin] = R.md_count[ibin]
Expand Down
2 changes: 1 addition & 1 deletion src/tools/contributions.jl
Original file line number Diff line number Diff line change
Expand Up @@ -83,7 +83,7 @@ function contributions(
warning = true,
)
(warning && s.nmols > 1) && warning_nmols_types()
indexes = [atom.index for atom in atoms]
indexes = PDBTools.index.(atoms)
# Check which types of atoms belong to this selection
selected_types = which_types(s, indexes, warning = warning)
return contributions(s, atom_contributions, selected_types)
Expand Down
38 changes: 28 additions & 10 deletions src/tools/coordination_number.jl
Original file line number Diff line number Diff line change
Expand Up @@ -73,6 +73,10 @@ function coordination_number(
return coordination_number(R, group_contributions)
end
function coordination_number(R::Result, group_contributions::Vector{Float64})
# obs: we accumulate group coordination numbers, but at the end we
# divide by random count of, so we the group contributions are the contributions
# to the MDDF. To recover the coordination number we need to multiply back
# the group contributions by the random count.
cn = cumsum(
group_contributions[i] * R.md_count_random[i] for
i in eachindex(group_contributions)
Expand All @@ -82,16 +86,30 @@ end

@testitem "coordination_number" begin
using ComplexMixtures, PDBTools
import ComplexMixtures.Testing: test_dir

pdb = readPDB("$test_dir/data/NAMD/structure.pdb")
R = load("$test_dir/data/NAMD/protein_tmao.json"; legacy_warning = false)
solute = Selection(PDBTools.select(pdb, "protein"), nmols = 1)
cn = coordination_number(solute, R.solute_atom, R, PDBTools.select(pdb, "residue 50"))
@test cn[findlast(<(5), R.d)] ≈ 0.24999999999999997 atol = 1e-10

pdb = readPDB("$test_dir/data/NAMD/Protein_in_Glycerol/system.pdb")
R = load("$test_dir/data/NAMD/Protein_in_Glycerol/protein_water.json"; legacy_warning=false)
using ComplexMixtures.Testing
dir = "$(Testing.data_dir)/NAMD"
atoms = readPDB("$dir/structure.pdb")
protein = Selection(select(atoms, "protein"), nmols = 1)
tmao = Selection(select(atoms, "resname TMAO"), natomspermol = 14)
options = Options(lastframe = 1, seed = 321, StableRNG = true, nthreads = 1, silent = true, n_random_samples=200)
traj = Trajectory("$dir/trajectory.dcd", protein, tmao)
R = mddf(traj, options)

# Test self-consistency
@test sum(R.solute_atom) ≈ sum(R.solvent_atom)
@test coordination_number(R, contributions(protein, R.solute_atom, select(atoms, "protein"))) == R.coordination_number

# Checked with vmd: same residue as (resname TMAO and within 3.0 of protein)
@test R.coordination_number[findfirst(>(3), R.d)] == 7.0
@test R.coordination_number[findfirst(>(5), R.d)] == 14.0

# THR93 forms a hydrogen bond with TMAO in this frame
thr93 = select(atoms, "protein and residue 93")
@test last(coordination_number(R, contributions(protein, R.solute_atom, thr93))) == 1

# Consistency of a read-out result
pdb = readPDB("$(Testing.data_dir)/NAMD/Protein_in_Glycerol/system.pdb")
R = load("$(Testing.data_dir)/NAMD/Protein_in_Glycerol/protein_water.json"; legacy_warning=false)
group = PDBTools.select(pdb, "protein")
solute = Selection(group, nmols = 1)
cn = coordination_number(solute, R.solute_atom, R, group)
Expand Down
2 changes: 1 addition & 1 deletion src/updatecounters.jl
Original file line number Diff line number Diff line change
Expand Up @@ -43,7 +43,7 @@ function updatecounters!(R::Result, system::AbstractPeriodicSystem, frame_weight
R.md_count_random[ibin] += frame_weight
if md.ref_atom_within_cutoff
ibin = setbin(md.d_ref_atom, R.options.binstep)
R.rdf_count_random[ibin] += 1
R.rdf_count_random[ibin] += frame_weight
end
end
return R
Expand Down
2 changes: 1 addition & 1 deletion test/data/Gromacs/EMI_DCA.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/Gromacs/EMI_EMI.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/Gromacs/protein_EMI.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/NAMD/protein_tmao.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/NAMD/tmao_tmao.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/NAMD/water_tmao.json

Large diffs are not rendered by default.

2 changes: 1 addition & 1 deletion test/data/NAMD/water_water.json

Large diffs are not rendered by default.