Bio.PDB.StructureBuilder module
Consumer class that builds a Structure object.
This is used by the PDBParser and MMCIFparser classes.
- class Bio.PDB.StructureBuilder.StructureBuilder
Bases:
object
Deals with constructing the Structure object.
The StructureBuilder class is used by the PDBParser classes to translate a file to a Structure object.
- __init__()
Initialize this instance.
- set_header(header)
Set header.
- set_line_counter(line_counter: int)
Tracks line in the PDB file that is being parsed.
- Arguments:
line_counter - int
- init_structure(structure_id: str)
Initialize a new Structure object with given id.
- Arguments:
structure_id - string
- init_model(model_id: int, serial_num: int | None = None)
Create a new Model object with given id.
- Arguments:
id - int
serial_num - int
- init_chain(chain_id: str)
Create a new Chain object with given id.
- Arguments:
chain_id - string
- init_seg(segid: str)
Flag a change in segid.
- Arguments:
segid - string
- init_residue(resname: str, field: str, resseq: int, icode: str)
Create a new Residue object.
- Arguments:
resname - string, e.g. “ASN”
field - hetero flag, “W” for waters, “H” for hetero residues, otherwise blank.
resseq - int, sequence identifier
icode - string, insertion code
- init_atom(name: str, coord: ndarray, b_factor: float, occupancy: float, altloc: str, fullname: str, serial_number=None, element: str | None = None, pqr_charge: float | None = None, radius: float | None = None, is_pqr: bool = False)
Create a new Atom object.
- Arguments:
name - string, atom name, e.g. CA, spaces should be stripped
coord - NumPy array (Float0, length 3), atomic coordinates
b_factor - float, B factor
occupancy - float
altloc - string, alternative location specifier
fullname - string, atom name including spaces, e.g. “ CA “
element - string, upper case, e.g. “HG” for mercury
pqr_charge - float, atom charge (PQR format)
radius - float, atom radius (PQR format)
is_pqr - boolean, flag to specify if a .pqr file is being parsed
- set_anisou(anisou_array)
Set anisotropic B factor of current Atom.
- set_siguij(siguij_array)
Set standard deviation of anisotropic B factor of current Atom.
- set_sigatm(sigatm_array)
Set standard deviation of atom position of current Atom.
- get_structure()
Return the structure.
- set_symmetry(spacegroup, cell)
Set symmetry.