Building DNA Structures¶
This guide covers all the ways to create DNA structures with MDNA: from simple sequences to custom-shaped assemblies.
Basic Construction¶
From a Sequence¶
From a Number of Base Pairs¶
When no sequence is provided, a random sequence is generated:
From Existing Data¶
Load from an MDTraj trajectory or PDB file:
import mdtraj as md
# From a PDB file
dna = mdna.load(filename='my_dna.pdb')
# From an MDTraj trajectory
traj = md.load('trajectory.xtc', top='topology.pdb')
dna = mdna.load(traj=traj, chainids=[0, 1])
Or from precomputed reference frames:
import numpy as np
frames = np.load('frames.npy') # shape (n_bp, 4, 3)
dna = mdna.load(frames=frames, sequence='ATCG...')
Custom Shapes¶
MDNA uses spline interpolation through control points to define arbitrary DNA paths.
Built-in Shapes¶
The Shapes class provides common geometries:
# Straight line (default)
points = mdna.Shapes.line(length=5)
# Circle
points = mdna.Shapes.circle(radius=3)
# Helix / supercoil
points = mdna.Shapes.helix(height=3, pitch=5, radius=7, num_turns=4)
Use them with make():
Custom Control Points¶
Define any 3D path with at least 4 points:
import numpy as np
control_points = np.array([
[0, 0, 0],
[30, 10, -10],
[50, 10, 20],
[3, 4, 30]
])
dna = mdna.make(n_bp=100, control_points=control_points)
dna.draw()
Tip
When only control_points are provided (no sequence or n_bp), MDNA infers n_bp from the spline arc length using the standard 0.34 nm rise per base pair.
Circular DNA¶
Create closed minicircles with optional linking number control:
# Relaxed minicircle
dna = mdna.make(n_bp=200, circular=True)
print('Lk, Wr, Tw:', dna.get_linking_number())
dna.draw()
# Supercoiled minicircle (overwound by 8 turns)
dna = mdna.make(n_bp=200, circular=True, dLk=8)
print('Lk, Wr, Tw:', dna.get_linking_number())
Note
When minimizing supercoiled DNA, use equilibrate_writhe=True to allow the writhe to equilibrate:
Minimization¶
Monte Carlo relaxation resolves steric clashes and elastic strain:
dna = mdna.make(n_bp=100)
dna.minimize(
temperature=300, # Kelvin
exvol_rad=2.0, # Excluded volume radius
endpoints_fixed=True # Pin the ends
)
You can also fix specific base pairs:
To view the Monte Carlo trajectory:
Extending DNA¶
Add base pairs to the 3' or 5' end of an existing structure:
dna = mdna.make(n_bp=50)
# Extend forward (3' direction)
dna.extend(sequence='G' * 30, forward=True)
# Extend backward (5' direction)
dna.extend(sequence='C' * 20, forward=False)
dna.draw()
Connecting Two DNA Structures¶
Join two Nucleic objects with an automatically optimized linker:
import numpy as np
dna0 = mdna.make(sequence='AAAAAAAAA')
dna1 = mdna.make(
sequence='GGGGGGGGG',
control_points=mdna.Shapes.line(1) + np.array([4, 0, -5])
)
# connect() finds the optimal number of linker base pairs
connected = mdna.connect(dna0, dna1)
connected.describe()
connected.draw()
The connect() function:
- Computes the twist difference between the two end frames
- Finds the optimal number of linking base pairs to achieve neutral twist
- Interpolates a connecting path between the two structures
- Optionally minimizes the resulting assembly
Exporting Structures¶
Save as PDB¶
Get MDTraj Trajectory¶
Direct PDB Generation¶
For one-shot PDB generation without creating a Nucleic object:
Summary¶
| Task | Function / Method |
|---|---|
| Generate from sequence/shape | mdna.make() |
| Load from trajectory/frames | mdna.load() |
| Relax structure | nucleic.minimize() |
| Extend DNA | nucleic.extend() |
| Join two structures | mdna.connect() |
| Export PDB | nucleic.save_pdb() |
| One-shot PDB | mdna.sequence_to_pdb() |
| MD simulation | mdna.sequence_to_md() |