Curvature Analysis
MembraneAnalysis.curvature
— Functioncurvature(;
point,
hq,
box_dims,
q_max,
q_min=0)
Calculates curvature per mode of a point from fluctuation spectrum.
Keyword arguments
point
: ordered pair of X and Y values;hq
: 2D matrix of height fluctuation spectrum;box_dims
: ordered pair of simulation box Lx and Ly values;q_max
: maximum q mode magnitude value to be used.q_min
: minimum q mode magnitude value to be used (zero by default).
MembraneAnalysis.lipids_sampled_curvature
— Functionlipids_sampled_curvature(;
pdb_file,
traj_files,
fs_files,
output_dir,
lipids,
q_max,
q_min=0)
Calculates mean sampled curvature of heavy atoms of each lipid. Results for lipid "XXXX" will be stored in XXXX_cs.dat
in the output directory.
Keyword arguments
pdb_file
: PDB structure file;traj_files
: a list of trajectory files;fs_files
: a list of corresponding HDF5 fluctuation spectrum files of the trajectory files;output_dir
: output directory;lipids
: a list of lipids of typeLipid
as defined inlipids.jl
;q_max
: maximum q mode magnitude value to be used.q_min
: minimum q mode magnitude value to be used (zero bu default).
MembraneAnalysis.peptide_sampled_curvature
— Functionpeptide_sampled_curvature(;
pdb_file,
traj_files,
fs_files,
output_dir,
lipids,
q_max)
Calculates mean sampled curvature of CA atoms of peptide residues.
Keyword arguments
pdb_file
: PDB structure file;traj_files
: a list of trajectory files;fs_files
: a list of corresponding HDF5 fluctuation spectrum files of the trajectory files;output_file
: output file;lipids
: a list of lipids of typeLipid
as defined inlipids.jl
;n_residues
: number of residues in the peptide;q_max
: maximum q mode magnitude value to be used.
MembraneAnalysis.lipids_curvature_spectrum
— Functionlipids_curvature_spectrum(;
pdb_file,
traj_files,
fs_files,
output_dir,
lipids,
ref_atoms=Dict(lipid => lipid.ref_atom for lipid in lipids),
q_max,
q_min,
nqs)
Calculates curvature spectrum of the lipids using their reference atom position. Assumes square (Lx = Ly) bilayer. Results for lipid "XXXX" will be stored in XXXX_cqs.dat
in the output directory.
Keyword arguments
pdb_file
: PDB structure file;traj_files
: a list of trajectory files;fs_files
: a list of corresponding HDF5 fluctuation spectrum files of the trajectory files;output_dir
: output directory;lipids
: a list of lipids of typeLipid
as defined inlipids.jl
;ref_atoms
: a dictionary of reference atoms for each lipid;q_max
: maximum q mode magnitude value to be used;q_min
: minimum q mode magnitude value to be used;nqs
: number of segments to divide the q range from qmin to qmax;tag_files
: a list of corresponding lipid interaction tag files of the trajectory files. Uses all lipids by default;tag_id
: the tag id of the selected lipids.
MembraneAnalysis.peptide_curvature_spectrum
— Functionpeptide_curvature_spectrum(;
pdb_file,
traj_files,
fs_files,
output_file,
lipids,
ref_residue),
q_max)
Calculates curvature spectrum of the peptide using the CA atom of its reference residue. Assumes square (Lx = Ly) bilayer.
Keyword arguments
pdb_file
: PDB structure file;traj_files
: a list of trajectory files;fs_files
: a list of corresponding HDF5 fluctuation spectrum files of the trajectory files;output_file
: output file;lipids
: a list of lipids of typeLipid
as defined inlipids.jl
;ref_residue
: residue number of the reference residue of the peptide;q_max
: maximum q mode magnitude value to be used.
MembraneAnalysis.TCB_analysis
— FunctionTCB_analysis(;
input_dir,
lipids,
weights=ones(length(lipids)) ./ length(lipids),
z_cutoff,
area=readdlm(input_dir * "A.dat")[1],
output_dir,
tcb_plot=false)
Calculates bilayer bending rigidity modulus and mean sampled curvature of lipids relative to a weighted average from transverse curvature bias analysis. Optionally plots mean sampled curvature of atoms of each lipid as a function of height.
Keyword arguments
input_dir
: directory with lipid atoms height and curvature files (e.g.XXXX_zs.dat
andXXXX_cs.dat
for lipid "XXXX");lipids
: a list of lipids of typeLipid
as defined inlipids.jl
;weights
: a list of the same size aslipids
, determining the weight of each lipid's TCB curve in the analysis. (Should be equal to the fraction of bilayer area covered by that lipid. Will be equal by default.);z_cutoff
: cutoff height to exclude anomalous behavior near lipid head region;area
: bilayer area, will be read fromA.dat
ininput_dir
by default;output_dir
: output directory;tcb_plot
: saves a plot (TCB_plot.pdf
) inoutput_dir
iftrue
.
MembraneAnalysis.hq2_analysis
— Functionhq2_analysis(;
input_dir,
n_points=size(readdlm(input_dir * "hq2.dat"))[1],
area=readdlm(input_dir * "A.dat")[1],
model="TD",
output_dir,
hq2_plot=false)
Calculates bilayer bending rigidity modulus by fitting |q|^4 × <|hq|^2> vs |q| data to the chosen model. Available models are "HC" (Helfrich-Canham), "MNK" (May-Narang-Kopelevich), or "TD" (Terzi-Deserno), as denoted by equations (36), (37), and (47) respectively in the following paper:
https://doi.org/10.1063/1.4990404
Saves the results to kc_hq2.dat
in output_dir
and optionally plots the fit.
Keyword arguments
input_dir
: input directory containinghq2.dat
;n_points
: number of |q| data points to use (uses all by default);model
: chosen theoretical model to fit <|hq|^2> vs |q|;area
: bilayer area, will be read fromA.dat
ininput_dir
by default;output_dir
: output directory;hq2_plot
: saves a plot (hq2_plot.pdf
) inoutput_dir
iftrue
.