next up previous contents index
Next: model.symmetry.define() define Up: Calculation of spatial restraints Previous: Specification of pseudo atoms   Contents   Index

model.restraints.make() -- make restraints

edat = <energy_data>   objective function parameters
aln = <alignment>   Template-model alignment; for homology-derived restraints only
io = <io_data>   Options for reading atom files
restraint_type = <str:1> 'STEREO' restraint type to be calculated: 'STEREO' | 'BOND' | 'ANGLE' | 'IMPROPER' | 'DIHEDRAL' | 'SPHERE' | 'SPHERE14' | 'LJ' | 'LJ14' | 'COULOMB' | 'COULOMB14' | 'ALPHA' | 'STRAND' | 'SHEET' | 'DISTANCE' | 'USER_DISTANCE' | 'NONB_PAIR_SPLINE' | 'PHI-PSI_BINORMAL' | 'PHI_DIHEDRAL' | 'PSI_DIHEDRAL' | 'OMEGA_DIHEDRAL' | 'CHI1_DIHEDRAL' | 'CHI2_DIHEDRAL' | 'CHI3_DIHEDRAL' | 'CHI4_DIHEDRAL'
dih_lib_only = <bool:1> False whether to use only library, not homologs for dihedral angle rsrs
mnch_lib = <int:1> 1 which MNCH lib to use in MAKE_RESTRAINTS
intersegment = <bool:1> True whether to restrain inter-segment non-bonded pairs
residue_grouping = <int:1> 1
maximal_distance = <float:1> 999 maximal distance for distance restraints
residue_span_range = <int:2> 0 99999 range of residues spanning the allowed distances; for MAKE_RESTRAINTS, PICK_RESTRAINTS, non-bonded dynamic pairs
residue_span_sign = <bool:1> True whether to do N*(N-1)/2 loop for atom pairs in MAKE_RESTRAINTS RESTRAINT_TYPE = 'distance'
restraint_sel_atoms = <int:1> 1 a restraint other than non-bonded pair has to have at least as many selected atoms
accessibility_type = <int:1> 8 type of solvent accessibility: 1-10
distance_rsr_model = <int:1> 1 the model for calculating distance restraints: 1-7
restraint_group = <int:1> 26 physical restraint group
restraint_stdev = <float:2> 0.0 1.0 transforming factors for standard deviations (y=a+bx) in models 1-6 or standard deviation for model 7 (a)
restraint_stdev2 = <float:3> 0 0 0 transforming standard deviation in models 3-6: S' = S + [ a + b max(0, c-g) ]
restraint_parameters = <float:0> 3 1 3 3 4 2 0 0.0 0.087 restraint parameters for 'USER_DISTANCE'
basis_pdf_weight = <str:1> 'LOCAL' a method for calculation of basis pdf weights: 'LOCAL' | 'GLOBAL'
basis_relative_weight = <float:1> 0.05 the cutoff weight of basis pdf's for their removal
residue_ids = <str:0> '' residue id (number:chnid)
spline_on_site = <bool:1> False whether to convert restraints to splines
spline_dx = <float:1> 0.5 interval size for splining restraints
spline_min_points = <int:1> 5 have at least as many intervals in a spline
spline_range = <float:1> 4.0 range of the splines
sheet_h_bonds = <int:1> 7 specify hydrogen bonds in a beta-sheet

Requirements:
topology & parameters [& picked atoms sets 2 and 3]

Description:
This command calculates and selects new restraints of a specified type. See the original papers for the most detailed definition and description of the restraints [Šali & Blundell, 1993,Šali & Overington, 1994]. The calculation of restraints of all types is now (partly) limited to the selected atoms only (either set 1, or 2 and 3; see below). The new restraints are added to any currently present.

restraint_type selects the types of the generated restraints. Only one restraint type can be selected at a time, except for the stereochemical restraints (BOND, ANGLE, DIHEDRAL, IMPROPER) that can all be calculated at the same time. It is useful to distinguish between the stereochemical restraints and homology-derived restraints. The stereochemical restraints are obtained from libraries that depend on atom and/or residue types only (e.g., CHARMM 22 force field [MacKerell et al., 1998] or statistical potentials), and do not require the alignment aln with template structures. In contrast, the homology-derived restraints are calculated from related protein structures, which correspond to all but the last sequence in the alignment aln (the target). These templates are read from coordinate files, which are the only data files required. All restraints are added to the existing restraints, even if they duplicate them (but see the comment for the 'OMEGA' restraints below).

Stereochemical restraints:

Homology-derived restraints:

For these restraints, the input alignment aln must be given.

basis_relative_weight is the cutoff for removing weak basis pdf's from poly-Gaussian feature pdf's: a basis pdf whose weight is less than the basis_relative_weight fraction of the largest weight is deleted.

Example: examples/commands/make_restraints.py


# Example for: model.restraints.make(), model.restraints.spline(),
#              model.restraints.write()

# This will compare energies of bond length restraints expressed
# by harmonic potential and by cubic spline.

log.verbose()
env = environ()
env.libs.topology.read(file='$(LIB)/top_heav.lib')
env.libs.parameters.read(file='$(LIB)/par.lib')

code = '1fas'
aln = alignment(env)
mdl = model(env, file=code, model_segment=('1:', '61:'))
aln.append_model(mdl, atom_files=code, align_codes=code)
aln.append_model(mdl, atom_files=code+'.ini', align_codes=code+'-ini')
mdl.generate_topology(aln, sequence=code+'-ini')
mdl.transfer_xyz(aln)
mdl.build(initialize_xyz=False, build_method='INTERNAL_COORDINATES')
mdl.write(file=code+'.ini')

mdl.restraints.make(aln, restraint_type='bond', spline_on_site=False)
mdl.restraints.write(file=code+'-1.rsr')
edat = energy_data(dynamic_sphere=False)
mdl.energy(edat=edat)

mdl.restraints.spline(spline_range=5.0, spline_dx=0.005,
                      spline_select=(3, 1, 1), edat=edat)
mdl.restraints.condense()
mdl.restraints.write(file=code+'-2.rsr')
mdl.energy(edat=edat)


next up previous contents index
Next: model.symmetry.define() define Up: Calculation of spatial restraints Previous: Specification of pseudo atoms   Contents   Index
Ben Webb 2006-02-28