openNucleome¶
openNucleome.whole_nucleus_model¶
-
class
OpenNucleome
(temperature=1.0, gamma=0.1, timestep=0.005, mass_scale=1.0)¶ The OpenNucleome class performes the whole nucleus dynamics based on the compartment annotations sequence of chromosomes.
The simulations can be performed using customed values for the type-to-type, intra, inter-chromosomal, chromosomes-nucleoli, and chromosomes-speckles parameters.
In our model, each component has specific type. Compartment A: 1, Compartment B: 2, Pericentromeres: 3, Centromeres: 4, Nucleoli: 5, dP Speckles: 6, P Speckles: 7, Lamina: 8
Initialize the whole nucleus model.
- temperature (float, required):
Temperature in reduced units. (Default value = 1.0).
- gamma (float, required):
Friction/Damping constant in units of reciprocal time (\(1/\tau\)) (Default value = 0.1).
- timestep (float, required):
Simulation time step in units of \(\tau\) (Default value = 0.005).
- mass_scale (float, required):
Mass scale used in units of \(\mu\) (Default value = 1.0).
-
add_chromosome_potential
(flag, ideal_file, types_file, inter_file)¶ Add the potential related to the genome in the system.
- flag (dict, required):
A dict contains the switch of each chromosome-related potential.
- ideal_file (string, required):
The path to the ideal potential scaling factor txt file.
- types_file (string, required):
The path to the type-to-type potential scaling factor txt file.
- inter_file (string, required):
The path to the inter-chromosomal potential scaling factor txt file.
-
add_lamina_potential
(flag, chr_lam_param)¶ Add the potential related to the nuclear lamina in the system.
- flag (dict, required):
A dict contains the switch of each lamina-related potential
- chr_lam_param (string, required):
The path to the potential scaling factor txt file.
-
add_nucleolus_potential
(flag, rescalar_file)¶ Add the potential related to the nucleoli in the system.
- flag (dict, required):
A dict contains the switch of each nucleolus-related potential
- rescalar_file (string, required):
The path to the txt file of rescaling factors (SPIN state probabilities).
-
add_speckle_potential
(flag, chr_spec_param)¶ Add the potential related to the nuclear speckles in the system.
- flag (dict, required):
A dict contains the switch of each speckle-related potential
- chr_spec_param (string, required):
The path to the potential scaling factor txt file.
-
construct_topology
(flag_membrane)¶ Construct the topology for the system.
This function is called automatically when creating a system.
- flag_membrane (bool, required):
A flag to control the dynamics of the lamina membrane, True: include the dynamics; False: exclude the dynamics (Default: False)
-
create_simulation
(platform_type='CPU')¶ Specify a platform when creating the simulation.
- platform_type (string, required):
The CPU platform is set as the default, while CUDA, OpenCL, and Reference are available as optional platforms.
-
create_system
(PDB_file, flag_membrane=False, lam_bond=None)¶ Create the system with initial configurations.
- PDB_file (string, required):
The path to the chromosome PDB file, including nucleoli, speckles, and lamina.
- flag_membrane (bool, required):
A flag to control the dynamics of the lamina membrane, True: include the dynamics; False: exclude the dynamics (Default: False)
- lam_bond (string, required):
The file containing the bond between specific lamina beads. (Default: None)
-
generate_element
(flag_membrane)¶ Generate elements for each coarse-grained bead.
This function is called automatically when creating a system.
- flag_membrane (bool, required):
A flag to control the dynamics of the lamina membrane, True: include the dynamics; False: exclude the dynamics (Default: False)
-
load_customized_settings
(force_field, param_folder, k=1.0)¶ Load customized force field settings
- force_field (Pandas Dataframe, required):
Store the flag of all the potentials and the parameter file names
- param_folder (string, required):
Save the customized parameters used in the simulations
- k (float, required):
The strength of force squeezing the nucleus (Default: 1.0)
-
load_default_settings
()¶ Load default force field settings
-
save_system
(system_xml)¶ Save the simulation state to readable xml format.
openNucleome.chromosome¶
-
class
Chromosome
(N_bead, N_chr, N_type, bond_list, compart_type, chr_groups, bead_groups, mol_type)¶ The Chromosome class stores all the useful potentials amoung chromatin beads used in the whole nucleus model.
Initialize the useful parameters
- N_bead (int, required):
the number of beads in the system.
- N_chr (int, required):
the number of chromosome beads in the system
- N_type (int, required):
the number of types in the system.
- bond_list (list, required):
contains the corresponding atom indexes of all bonds.
- compart_type (list, required):
contains the corresponding compartment types of all beads.
- chr_groups (list, required):
contains the indexes of all chromatin beads.
- bead_groups (list, required):
contains the indexes of beads of each type
- mol_type (list, required):
contains the indexes of all molecules (each chromosome is one molecule, but each neclear body bead is also one molecule).
-
add_angle_force
(ka=2.0)¶ Add an angular potential between bonds connecting beads.
- ka (float, required):
angle coefficients, default = 2.0
-
add_class2_bond
(force_group, r0c=0.5, kc2=20.0, kc3=20.0, kc4=20.0)¶ Internal function that inits Class2 bond force. The default parameters value is for the 100KB resolution. (Recommended)
- force_group (int, required):
labels the index of the current force.
- r0c, kc2, kc3, kc4 (float, required):
class2 bond coefficients
-
add_fene_bond
(force_group, kfb=30.0, r0=1.5, epsilon_fene=1.0, sigma_fene=1.0, cutoff_fene=1.12)¶ Internal function that inits FENE bond force. The default parameters value is for the 1MB resolution.
- force_group (int, required):
labels the index of the current force.
- r0c, kfb, r0, epsilon_fene, sigma_fene, cutoff_fene (float, required):
fene bond coefficients
-
add_ideal_potential
(ideal_file, force_group, rc_tanh_ideal=0.54, dend_ideal=1000, cutoff_ideal=2.0)¶ Add the Intra-chromosomal ideal potential using custom values for interactions between beads separated by a genomic distance \(d\).
The parameter rc_tanh_ideal is part of the probability of crosslink function
\[\begin{split}f(r_{i,j}) & = \frac{1}{2}\left(1+\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right) \\ & + \frac{1}{2}\left(1-\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right)\left(\frac{r_c}{r_{i,j}}\right)^4,\end{split}\]where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- ideal_file (string, required):
The path to the ideal potential scaling factor txt file.
- force_group (int, required):
labels the index of the current force.
- dend_ideal (int, required):
Cutoff of genomic separation (Default value = 1000 for the 100KB model).
- rc_tanh_ideal (float, required):
\(r_c\) in the equation, a good estimation of the bond length during the simulation (Default value = 0.54 for the 100KB model).
- cutoff_ideal (float, required):
Cutoff of Ideal potential (Default value = 2.0 for the 100KB model).
-
add_inter_potential
(inter_file, force_group, rc_tanh_inter=0.54, cutoff_inter=2.0)¶ Add the Inter-chromosome Chromosome potential using custom values for interactions between beads separated by a genomic distance \(d\).
The parameter rc_tanh_inter is part of the probability of crosslink function
\[\begin{split}f(r_{i,j}) & = \frac{1}{2}\left(1+\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right) \\ & + \frac{1}{2}\left(1-\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right)\left(\frac{r_c}{r_{i,j}}\right)^4,\end{split}\]where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- inter_file (string, required):
The path to the inter-chromosomal potential scaling factor txt file.
- force_group (int, required):
labels the index of the current force.
- rc_tanh_inter (float, required):
\(r_c\) in the equation, a good estimation of the bond length during the simulation (Default value = 0.54 for the 100KB model).
- cutoff_inter (float, required):
Cutoff of inter-chromosomal potential (Default value = 2.0 for the 100KB model).
-
add_softcore
(force_group, epsilon_sc=1.0, sigma_sc=0.5, cutoff_sc=0.56, Ecut_sc=4.0)¶ We plug the tanh potential into the WCA potential to represent the excluded volume effect;
The softcore potential will only be applied to the chromosomes (type 1, 2, 3, 4).
When \(r_{sc} < r < r_{cutoff}\), use \(E_{LJ} = 4\epsilon_{sc}\left(\left(\frac{\sigma_{sc}}{r}\right)^{12} -\left(\frac{\sigma_{sc}}{r}\right)^{6} \right) + \epsilon_{rc}\); \(r < r_{sc}\), use \(E = \frac{1}{2}E_{cut}\left(1+\text{tanh}(\frac{2E_{LJ}}{E_{cut}}-1) \right)\), where \(r_{sc} \approx 0.485\).
- force_group (int, required):
labels the index of the current force.
- epsilon_sc (float, required):
energy units (Default value = 1.0).
- sigma_sc (float, required):
distance units, which represents the diameter of a 100KB-resolution bead (Default value = 0.5).
- cutoff_sc (float, required):
distance units, the cutoff (Default value = 0.5*1.12).
- Ecut_sc (float, required):
energy units, the energy cost for the chain passing (Default value = 4.0).
-
add_type_type_potential
(types_file, force_group, rc_tanh_type=0.54, cutoff_type=2.0)¶ Add the type-to-type potential using custom values for interactions between the chromatin types.
The parameter rc_tanh_type is part of the probability of crosslink function
\[\begin{split}f(r_{i,j}) & = \frac{1}{2}\left(1+\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right) \\ & + \frac{1}{2}\left(1-\text{tanh}\left[(r_c - r_{i,j})^{-5}+5(r_c - r_{i,j})\right] \right)\left(\frac{r_c}{r_{i,j}}\right)^4,\end{split}\]where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- types_file (string, required):
The path to the type-to-type potential scaling factor txt file.
- force_group (int, required):
labels the index of the current force.
- rc_tanh_type (float, required):
\(r_c\) in the equation, a good estimation of the bond length during the simulation (Default value = 0.54 for the 100KB model).
- cutoff_type (float, required):
Cutoff of Type-Type potential (Default value = 2.0 for the 100KB model).
openNucleome.speckle¶
-
class
Speckle
(N_bead, N_type, bond_list, compart_type, chr_groups, bead_groups, mol_type)¶ The Speckle class stores all the useful potentials related to speckles used in the whole nucleus model.
Initialize the useful parameters
- N_bead (int, required):
the number of beads in the system.
- N_type (int, required):
the number of types in the system.
- bond_list (list, required):
contains the corresponding atom indexes of all bonds.
- compart_type (list, required):
contains the corresponding compartment types of all beads.
- chr_groups (list, required):
contains the indexes of all chromatin beads.
- bead_groups (list, required):
contains the indexes of beads of each type
- mol_type (list, required):
contains the indexes of all molecules (each chromosome is one molecule, but each neclear body bead is also one molecule).
-
add_chr_spec
(chr_spec_param, force_group, sigma_tanh_chr_spec=4.0, rc_tanh_chr_spec=0.75, cutoff_chr_spec=2.0)¶ Add nonbonded potential using custom values for interactions between chromosomes and speckles
The parameters are part of the probability of function \(g(r_{i,j}) = \frac{1}{2}\left(1 + \text{tanh}\left[\sigma(r_c - r_{i,j})\right] \right)\), where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- chr_spec_param (string, required):
The path to the potential scaling factor txt file.
- force_group (int, required):
labels the index of the current force.
- sigma_tanh_chr_spec (float, required):
Distance units (Default value = 4.0).
- rc_tanh_chr_spec (float, required):
Distance units, where \(g(r_{i,j}) = 0.5\) (Default value = 0.75).
- cutoff_chr_spec (float, required):
The cutoff (Default value = 2.0)
-
add_spec_spec
(force_group, epsilon_dP_dP=3.0, sigma_spec=0.5, cutoff_dP_dP=1.5, epsilon_P=1.0, cutoff_P=0.56)¶ Add nonbonded rescaled LJpotential between Speckles-Speckles
Notice that the system has two different speckles dP (type 6) and P (type 7).
- force_group (int, required):
labels the index of the current force.
- epsilon_dP_dP (float, required):
Epsilon of LJ between (dP,dP), energy units (Default value = 3.0).
- sigma_spec (float, required):
Sigma of LJ between speckles and speckles, distance units (Default value = 0.5).
- cutoff_dP_dP (float, required):
Cutoff of LJ between (dP,dP), distance units (Default value = 1.5).
- epsilon_P (float, required):
Epsilon of LJ between (dP,P) and (P,P), energy units (Default value = 1.0).
- cutoff_P (float, required):
Cutoff of LJ between (dP,P) and (P,P), distance units (Default value = 0.5*1.12).
openNucleome.nucleolus¶
-
class
Nucleolus
(N_bead, N_type, N_chr, N_chr_nuc, bond_list, compart_type, chr_groups, bead_groups, mol_type)¶ The Nucleolus class stores all the useful potentials related to nucleoli used in the whole nucleus model.
Initialize the useful parameters
- N_bead (int, required):
the number of beads in the system.
- N_type (int, required):
the number of types in the system.
- bond_list (list, required):
contains the corresponding atom indexes of all bonds.
- compart_type (list, required):
contains the corresponding compartment types of all beads.
- chr_groups (list, required):
contains the indexes of all chromatin beads.
- bead_groups (list, required):
contains the indexes of beads of each type
- mol_type (list, required):
contains the indexes of all molecules (each chromosome is one molecule, but each neclear body bead is also one molecule).
-
add_chr_nuc
(rescalar_file, force_group, sigma_tanh_chr_nuc=4.0, rc_tanh_chr_nuc=0.75, cutoff_chr_nuc=2.0)¶ Add nonbonded potential using custom values for interactions between chromosomes and nucleolus.
The parameters are part of the probability of function \(g(r_{i,j}) = \frac{1}{2}\left(1 + \text{tanh}\left[\sigma(r_c - r_{i,j})\right] \right)\), where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- rescalar_file (string, required):
The path to the txt file of rescaling factors (SPIN state probabilities).
- force_group (int, required):
labels the index of the current force.
- sigma_tanh_chr_nuc (float, required):
Distance units (Default value = 4.0).
- rc_tanh_chr_nuc (float, required):
Distance units, where \(g(r_{i,j}) = 0.5\) (Default value = 0.75).
- cutoff_chr_nuc (float, required):
The cutoff (Default value = 2.0).
-
add_nuc_nuc
(force_group, epsilon_nuc=3.0, sigma_nuc=0.5, cutoff_nuc=1.5)¶ Add nonbonded LJpotential between Nucleoli-Nucleoli.
- force_group (int, required):
labels the index of the current force.
- epsilon_nuc (float, required):
energy units (Default value = 3.0).
- sigma_nuc (float, required):
distance units (Default value = 0.5).
- cutoff_nuc (float, required):
The cutoff, distance units (Default value = 1.5).
-
add_nuc_spec
(force_group, epsilon_LJ_plain=1.0, sigma_LJ_plain=0.5, cutoff_LJ_plain=0.56)¶ Add nonbonded plain LJpotential between Nucleoli and Speckles
- force_group (int, required):
labels the index of the current force.
- epsilon_LJ_plain (float, required):
energy units (Default value = 1.0).
- sigma_LJ_plain (float, required):
distance units (Default value = 0.5).
- cutoff_LJ_plain (float, required):
The cutoff, distance units (Default value = 0.5*1.12).
openNucleome.lamina¶
-
class
Lamina
(N_bead, N_type, N_chr_nuc_spec, bond_list, compart_type, chr_groups, bead_groups, mol_type)¶ The Lamina class stores all the useful potentials related to lamina used in the whole nucleus model.
Initialize the useful parameters
- N_bead (int, required):
the number of beads in the system.
- N_type (int, required):
the number of types in the system.
- bond_list (list, required):
contains the corresponding atom indexes of all bonds.
- compart_type (list, required):
contains the corresponding compartment types of all beads.
- chr_groups (list, required):
contains the indexes of all chromatin beads.
- bead_groups (list, required):
contains the indexes of beads of each type
- mol_type (list, required):
contains the indexes of all molecules (each chromosome is one molecule, but each neclear body bead is also one molecule).
-
add_chr_lam
(chr_lam_param, force_group, sigma_tanh_chr_lam=4.0, rc_tanh_chr_lam=0.75, cutoff_chr_lam=2.0)¶ Add nonbonded potential using custom values for interactions between chromosomes and lamina
The parameters are part of the probability of function \(g(r_{i,j}) = \frac{1}{2}\left(1 + \text{tanh}\left[\sigma(r_c - r_{i,j})\right] \right)\), where \(r_{i,j}\) is the spatial distance between loci (beads) i and j.
- chr_lam_param (string, required):
The path to the potential scaling factor txt file.
- force_group (int, required):
labels the index of the current force.
- sigma_tanh_chr_lam (float, required):
Distance units (Default value = 4.0).
- rc_tanh_chr_lam (float, required):
Distance units, where \(g(r_{i,j}) = 0.5\) (Default value = 0.75).
- cutoff_chr_lam (float, required):
The cutoff (Default value = 2.0)
-
add_hardwall
(force_group, epsilon_hw=1.0, sigma_hw=0.75, cutoff_hw=0.84)¶ Add nonbonded hard-wall potential between (chromosomes, speckles, nucleoli) and lamina
The size of lamina is 1.0 (diameter) and the size of other beads is 0.5, so the sigma is 1.0/2 + 0.5/2 = 0.75 for the 100KB model
- force_group (int, required):
labels the index of the current force.
- epsilon_hw (float, required):
energy units (Default value = 1.0).
- sigma_hw (float, required):
distance units (Default value = 0.75).
- cutoff_hw (float, required):
The cutoff, distance units (Default value = 0.75*1.12).
-
add_lam_lam
(force_group, epsilon_lam=1.0, sigma_lam=0.5, cutoff_lam=0.5612310241546865)¶ Add nonbonded plain LJpotential between Lamina and Lamina
- force_group (int, required):
labels the index of the current force.
- epsilon_lam (float, required):
energy units (Default value = 1.0).
- sigma_lam (float, required):
Distance units (Default value = 0.5).
- cutoff_lam (float, required):
The cutoff (Default value = 0.5*1.12)
-
add_squeeze_nucleus
(force_group, k=1.0, R=13.0)¶ Add external force to squeeze the cell nucleus (Default: spring force)
- force_group (int, required):
labels the index of the current force.
- k (float, required):
energy units (Default value = 1.0).
- R (float, required):
perfect sphere nucleus radius (Default value = 13.0)
openNucleome.utils.DamID_TSASeq_calculation¶
-
DamID_TSASeq_calculation
(traj_data, gLength, maternalIdx, paternalIdx, start_frame, end_frame)¶ The function computes the DamID and TSA-Seq signals. The returns of this function are DamID, TSA-Seq, and number of speckle clusters.
- traj_data (string, required):
The trajectory file
- gLength (string, required):
The difference between each two neighboring values in the file represents the length of the chromosome
- maternalIdx (string, required):
Index of each maternal chromosome
- paternalIdx (string, required):
Index of each paternal chromosome
- start_frame (int, required):
The starting analyzed frame
- end_frame (int, required):
The ending analyzed frame
openNucleome.utils.coor_transformation¶
-
coor_transformation
(traj_data, output_traj)¶ The function converts the trajectory with position in OpenMM default unit to reduced unit
- traj_data (string, required):
The trajectory file from a simulation
- output_traj (string, required):
The trajectory with the reduced units
openNucleome.utils.final_frame¶
-
final_frame
(traj_data, type_final, output_config, n_chrom=46, n_nuc=300, n_spec=1600, nuc_type=5, lam_type=8)¶ The function extracts the final frame of the previous simulation to restart a new simulation.
- traj_data (string, required):
The trajectory file from a simulation
- type_final (string, required):
The file contains the type of each bead in the last frame
- output_config (string, required):
The final frame of the trajectory
- n_chrom (int, required):
The number of chromosomes (Default value = 46)
- n_nuc (int, required):
The number of nucleolus (Default value = 300)
- n_spec (int, required):
The number of speckles (Default value = 1600)
- nuc_type (int, required):
The type index of nucleolus (Default value = 5)
- lam_type (int, required):
The type index of lamina (Default value = 8)