io - Input and Output of Sequences, Structures and Maps

The io module deals with the input and output of entities, alignments, sequences, images. Importers for common file formats containing molecules such as PDB, SDF and CHARMM trajectory files are available. Sequence and alignment file formats such as FASTA and CLUSTALW are supported as well as various image data (e.g. png, dm3) and density map files (e.g. CCP4, MRC).

Molecular Structures

Loading Molecular Structures

The io modules offer several ways to load molecular structures depending on your requirements. The most general way is offered by LoadEntity(), which will automatically detect the file format based on the file extension.

LoadEntity(filename, format='auto')

Load entity from disk. If format is set to ‘auto’, the function guesses the filetype based on the extension of the file e.g. files ending in ‘.pdb’, ‘.ent’, ‘.ent.gz’ and ‘.pdb.gz’ will automatically be loaded as PDB files. For files without or exotic extensions, the format can be set explicitly as the second parameter.

# Recognizes SDF file by file extension
ent = io.LoadEntity('file.sdf')

# In this case, there is no file extensions, so you have to say it's a
# SDF file explicitly
ent = io.LoadEntity('file', 'sdf')

For a list of file formats supported by LoadEntity(), see Supported Structure File Formats.

Raises:

IOUnknownFormatException if the format string supplied is not recognized or the file format can not be detected based on the file extension.

IOException if the import fails due to an erroneous or inexistent file.

Some of the formats have a dedicated function that allows you to tweak many parameters that affect the import. PDB files can be loaded with LoadPDB() and mmCIF files with LoadMMCIF() (this also gives you access to the MMCifInfo class). It offers tighter control over the exact loading behaviour.

LoadPDB(filename, restrict_chains='', no_hetatms=None, fault_tolerant=None, load_multi=False, join_spread_atom_records=None, calpha_only=None, profile='DEFAULT', remote=False, remote_repo='pdb', dialect=None, seqres=False, bond_feasibility_check=None, read_conect=False)

Load PDB file from disk and return one or more entities. Several options allow to customize the exact behaviour of the PDB import. For more information on these options, see IO Profiles for entity importer.

Residues are flagged as ligand if they are mentioned in a HET record.

Parameters:
  • filename (str) – File to be loaded

  • restrict_chains (str) – If not an empty string, only chains listed in the string will be imported.

  • no_hetatms (bool) – If set to True, HETATM records will be ignored. Overrides the value of IOProfile.no_hetatms

  • fault_tolerant (bool) – Enable/disable fault-tolerant import. If set, overrides the value of IOProfile.fault_tolerant.

  • load_multi (bool) – If set to True, a list of entities will be returned instead of only the first. This is useful when dealing with multi-PDB files.

  • join_spread_atom_records (bool) – If set, overrides the value of IOProfile.join_spread_atom_records.

  • calpha_only (bool) – When set to true, forces the importer to only load atoms named CA. If set, overrides the value of IOProfile.calpha_only.

  • profile (str/ost.io.IOProfile) – Aggregation of flags and algorithms to control import and processing of molecular structures. Can either be a str specifying one of the default profiles [‘DEFAULT’, ‘SLOPPY’, ‘CHARMM’, ‘STRICT’] or an actual object of type ost.io.IOProfile. If a str defines a default profile, IOProfile.processor is set to ost.conop.RuleBasedProcessor with the currently set ost.conop.CompoundLib available as ost.conop.GetDefaultLib(). If no ost.conop.CompoundLib is available, ost.conop.HeuristicProcessor is used instead. See IO Profiles for entity importer for more info.

  • remote (bool) – If set to True, the method tries to load the pdb from the remote repository given as remote_repo. The filename is then interpreted as the entry id as further specified for the remote_repo parameter.

  • remote_repo (str) – Remote repository to fetch structure if remote is True. Must be one in [‘pdb’, ‘smtl’, ‘pdb_redo’]. In case of ‘pdb’ and ‘pdb_redo’, the entry must be given as lower case pdb id, which loads the deposited assymetric unit (e.g. ‘1ake’). In case of ‘smtl’, the entry must also specify the desired biounit (e.g. ‘1ake.1’).

  • dialect (str) – Specifies the particular dialect to use. If set, overrides the value of IOProfile.dialect

  • seqres (bool) – Whether to read SEQRES records. If set to True, the loaded entity and seqres entry will be returned as a tuple. If file doesnt contain SEQRES records, the returned ost.seq.SequenceList will be invalid.

  • bond_feasibility_check (bool) – Flag for IOProfile.processor. If turned on, atoms are only connected by bonds if they are within a reasonable distance (as defined by ost.conop.IsBondFeasible()). If set, overrides the value of ost.conop.Processor.check_bond_feasibility

  • read_conect (bool) – By default, OpenStructure doesn’t read CONECT statements in a pdb file. Reason is that they’re often missing and we prefer to rely on the chemical component dictionary from the PDB. However, there may be cases where you really want these CONECT statements. For example novel compounds with no entry in the chemical component dictionary. Setting this to True has two effects: 1) CONECT statements are read and blindly applied 2) The processor does not connect any pair of HETATOM atoms in order to not interfer with the CONECT statements.

Return type:

EntityHandle or a list thereof if load_multi is True.

Raises:

IOException if the import fails due to an erroneous or inexistent file

LoadMMCIF(filename, fault_tolerant=None, calpha_only=None, profile='DEFAULT', remote=False, seqres=False, info=False)

Load a mmCIF file and return one or more entities. Several options allow to customize the exact behaviour of the mmCIF import. For more information on these options, see IO Profiles for entity importer.

Residues are flagged as ligand if they are not waters nor covered by an entity_poly record (ie. they are non-polymer entities in pdbx_entity_nonpoly). Note that all residues except waters will be flagged as ligands if seqres=False (the default).

Parameters:
  • filename (str) – File to be loaded

  • fault_tolerant (bool) – Enable/disable fault-tolerant import. If set, overrides the value of IOProfile.fault_tolerant.

  • calpha_only (bool) – When set to true, forces the importer to only load atoms named CA. If set, overrides the value of IOProfile.calpha_only.

  • profile (str/ost.io.IOProfile) – Aggregation of flags and algorithms to control import and processing of molecular structures. Can either be a str specifying one of the default profiles [‘DEFAULT’, ‘SLOPPY’, ‘CHARMM’, ‘STRICT’] or an actual object of type ost.io.IOProfile. If a str defines a default profile, IOProfile.processor is set to ost.conop.RuleBasedProcessor with the currently set ost.conop.CompoundLib available as ost.conop.GetDefaultLib(). If no ost.conop.CompoundLib is available, ost.conop.HeuristicProcessor is used instead. See IO Profiles for entity importer for more info.

  • remote (bool) – If set to True, the method tries to load the pdb from the remote pdb repository www.pdb.org. The filename is then interpreted as the pdb id.

  • seqres (bool) – Whether to return SEQRES records. If True, a SequenceList object is returned as the second item. The sequences in the list are named according to the mmCIF chain name. This feature requires a default compound library to be defined and accessible via GetDefaultLib(). One letter codes of non standard compounds are set to X otherwise.

  • info (bool) – Whether to return an info container with the other output. If True, a MMCifInfo object is returned as last item.

Return type:

EntityHandle (or tuple if seqres or info are True).

Raises:

IOException if the import fails due to an erroneous or non-existent file.

PDBStrToEntity(pdb_string, profile=IOProfile(), process=False)

Load entity from a string in PDB format. By default the entity is loaded with an empty IO Profile and is not processed with the Processor, even if one is available in the IO Profile.

Parameters:
  • pdb_string – A PDB file as a string.

  • profile – The IO Profile to read the entity with. For more information on the IO Profiles available, see IO Profiles for entity importer.

  • process – If set to True, run the Processor contained in the IO Profile.

Return type:

EntityHandle.

To get an entity equivalent to one loaded with LoadPDB(), set the profile and process arguments as follows:

with open('protein.pdb') as pdb_fd:
    pdb_str = pdb.read()
    ent = io.PDBStrToEntity(pdb_str, ost.io.profiles['DEFAULT'], True)
LoadSDF(filename)

Load an SDF file and return an entity.

Parameters:

filename (str) – File to be loaded

Return type:

EntityHandle

SDFStrToEntity(sdf_string)

Load entity from a string in SDF format.

Parameters:

pdb_string – A SDF file as a string.

Return type:

EntityHandle.

class ost.io.OMF

Experimental data format capable of storing minimally required information of a ost.mol.EntityHandle in a compressed manner (OMF - OpenStructure Minimal Format). Shares lots of ideas with the mmtf or binaryCIF formats but comes with no dependencies attached. User can provide an upper error bound of atom coordinates in A which has an effect on compression ratio.

static FromEntity(ent, max_error=0.0, options=0)

Generates ost.io.OMF object starting from an ost.mol.EntityHandle.

Parameters:
  • ent (ost.mol.EntityHandle) – Structure data

  • max_error (float) – Maximum error in A of stored atom coordinates - affects compression ratio

  • options (ost.io.OMFOption) – Control file compression

Returns:

The created ost.io.OMF object

static FromFile(filepath)

Generates ost.io.OMF object from a file stored with ToFile() or a byte string created from ToBytes() that has been dumped to disk.

Parameters:

filepath (str) – The file

Returns:

The created ost.io.OMF object

static FromBytes(bytes_string)

Generates ost.io.OMF object from a bytes object created with ToBytes().

Parameters:

bytes_string (bytes) – The super heavily compressed OMF structure

ToFile(filepath)

Stores object to file

Parameters:

filepath (str) – The file

ToBytes()

Stores object to bytes string

Returns:

The bytes string

GetMaxError()

Get maximum error parameter

Returns:

The max error in A of the stored coordinates

GetAU()

Getter for assymetric unit

Returns:

The assymetric unit as ost.mol.EntityHandle

GetEntity()

Same as GetAU()

GetAUChain(chain_name)

Getter for single chain in the assymetric unit

Parameters:

chain_name (str) – The name of the chain you want

Returns:

ost.mol.EntityHandle of the AU only containing the specified chain

GetEntityChain(chain_name)

Same as GetAUChain()

GetName()
Returns:

The entity name, accessible as ost.mol.EntityHandle.GetName from the entity that has been used to create OMF object

GetChainNames()
Returns:

list with chain names of the stored structure

class ost.io.OMFOption

OMF options control file compression. They can be combined with bitwise or:

omf_opts = io.OMFOption.DEFAULT_PEPLIB | io.OMFOption.AVG_BFACTORS
  • DEFAULT_PEPLIB: Compound information (connectivity, chem type etc.) for each unique residue in a structure is stored in a library. This option defines a default library of such compounds with proteinogenic amino acids. This reduces the size of the library that needs to be written for each OMF file. IMPORTANT: This option has an impact on compression ratio as internal coordinate parameters are only available through the default peptide library

  • AVG_BFACTORS: Protein models often have the same bfactor for all atoms of a residue. One example is AFDB. Enabling this option results in only one bfactor being stored per residue, thus reducing file size.

  • ROUND_BFACTORS: Round bfactor to integer. This gives sufficient accuracy for bfactors in AFDB that are in range [0.0, 100.0].

  • SKIP_SS: Dont save secondary structure to save some memory

  • INFER_PEP_BONDS: Connectivity of compounds is stored in an internal library. However, inter-residue connectivity needs to be stored separately. If this option is enabled, stereochemically reasonable peptide bonds are inferred at loading and there is no need to save them.

Loading Molecular Structures From Remote Repositories

LoadPDB() already provides access to selected structural databases in the internet when enabling the remote flag. Predefined ost.io.remote.RemoteRepository objects are available as

from ost.io import remote
repo_name = 'smtl'
repo = remote.REMOTE_REPOSITORIES.get(repo_name)

# url for entry with id 1ake.1
print(repo.URLForID('1ake.1'))

where repo_name can be one of [‘pdb’, ‘cif’, ‘pdb_redo’, ‘smtl’]. Instead of explicit access, you can directly fetch data using:

RemoteGet(id, from_repo='pdb')

Invokes RemoteRepository.Get() on predefined repositories (‘pdb’, ‘smtl’, ‘cif’, ‘pdb_redo’)

Parameters:

from_repo (str) – One of the predefined repositories

RemoteLoad(id, from_repo='pdb')

Invokes RemoteRepository.Load() on predefined repositories (‘pdb’, ‘smtl’, ‘cif’, ‘pdb_redo’)

Parameters:

from_repo (str) – One of the predefined repositories

class RemoteRepository(name, url_pattern, type, id_transform='upper')

A remote repository represents a structural database accessible through the internet, e.g. the PDB or SWISS-MODEL template library.

Parameters:
  • name (str) – Name of the repository

  • url_pattern (str) – URL pattern for repository. Required format is described in URLForID()

  • type (str) – Data format to expect at resolved URL must be in (‘pdb’, ‘cif’)

  • id_transform (str) – Transformation to apply to ID before resolving URL in URLForID(). Must be in (‘lower’, ‘upper’)

Get(id)

Resolves URL with URLForID(), dumps the content in a temporary file and returns its path.

Parameters:

id (str) – ID to resolve

Load(id)

Resolves URL with URLForID() and directly loads/returns the according ost.mol.EntityHandle. Loading invokes the ost.io.LoadPDB()/ost.io.LoadMMCIF() with default parameterization. If you need custom settings, you might want to consider to call Get() and do the loading manually.

Parameters:

id (str) – ID to resolve

URLForID(id)

Resolves URL given url_pattern and id_transform provided at object initialization. The url_pattern must contain substring ‘$ID’. Given id, the URL to the structure gets constructed by applying id_transform and inserting it at the location of ‘$ID’. e.g. ‘https://files.rcsb.org/view/$ID.pdb’ given 1ake as id and ‘upper’ as id_transform resolves to: ‘https://files.rcsb.org/view/1AKE.pdb

Saving Molecular Structures

Saving a complete entity or a view is a matter of calling SaveEntity().

ent = io.LoadEntity('protein.pdb')
# save full entity
io.SaveEntity(ent, 'full.pdb')
# only save C-alpha atoms
io.SaveEntity(ent.Select('aname=CA and peptide=true'), 'calpha.pdb')

SavePDB() provides a simple way to save several entities into one file:

ent = io.LoadEntity('protein.pdb')
# Save complete entity
io.SavePDB(ent, 'full.pdb')
# Save chain A and chain B separately
io.SavePDB([ent.Select('cname=A'), ent.Select('cname=B')], 'split.pdb')
SaveEntity(ent, filename, format='auto')

Save entity to disk. If format is set to ‘auto’, the function guesses the filetype based on the file extension, otherwise the supplied format is checked against the available export plugins.

Parameters:
  • ent (EntityHandle or EntityView) – The entity to be saved

  • filename (string) – The filename

  • format (string) – Name of the format

Raises:

IOUnknownFormatException if the format string supplied is not recognized or the file format can not be detected based on the file extension.

SavePDB(models, filename, dialect=None, pqr=False, profile='DEFAULT')

Save entity or list of entities to disk. If a list of entities is supplied the PDB file will be saved as a multi PDB file. Each of the entities is wrapped into a MODEL/ENDMDL pair.

If the atom number exceeds 99999, ‘*’ is used.

Parameters:
  • models – The entity or list of entities (handles or views) to be saved

  • filename (string) – The filename

Raises:

IOException if the restrictions of the PDB format are not satisfied (with the exception of atom numbers, see above):

  • Chain names with more than one character

  • Atom positions with coordinates outside range [-999.99, 9999.99]

  • Residue names longer than three characters

  • Atom names longer than four characters

  • Numeric part of ost.mol.ResNum outside range [-999, 9999]

  • Alternative atom indicators longer than one character

SaveMMCIF(ent, filename, compound_lib=<ost.conop._ost_conop.CompoundLib object>, data_name='OST_structure', mmcif_conform=True, entity_info=<ost.io._ost_io.MMCifWriterEntityList object>)

Save OpenStructure entity in mmCIF format

Parameters:
  • ent (ost.mol.EntityHandle/ost.mol.EntityView) – OpenStructure Entity to be saved

  • filename (str) – Filename - .gz suffix triggers gzip compression

  • compound_lib (ost.conop.CompoundLib) – Compound library required when writing, uses ost.conop.GetDefaultLib() if not given

  • data_name (str) – Name of data block that will be written to mmCIF file. Typically, thats the PDB ID or some identifier.

  • mmcif_conform (bool) – Controls processing of structure, i.e. identification of mmCIF entities etc. before writing. Detailed description in mmCIF Writer. In short: If mmcif_conform is set to True, Chains in ent are expected to be valid mmCIF entities with residue numbers set according to underlying SEQRES. That should be the case when ent has been loaded with LoadMMCIF(). If mmcif_conform is set to False, heuristics kick in to identify and separate mmCIF entities based on ost.mol.ChemClass of the residues in a chain.

  • entity_info (MMCifWriterEntityList) – Advanced usage - description in mmCIF Writer

EntityToPDBStr(ent, profile=IOProfile())

Return entity as a string in PDB format.

Parameters:
Raises:

IOException if the restrictions of the PDB format are not satisfied (see ost.io.SavePDB())

Return type:

string.

EntityToSDFStr(ent)

Return entity as a string in SDF format.

Parameters:

entity – The EntityHandle or EntityView

Return type:

string.

SaveSDF(ent, filename)

Save entity to disk as an SDF file.

Parameters:

Sequences and Alignments

Loading sequence or alignment files

LoadSequence(filename, format='auto')

Load sequence data from disk. If format is set to ‘auto’, the function guesses the filetype based on the extension of the file. Files ending in ‘.fasta’, ‘.aln’ will automatically be loaded.

Parameters:
  • filename (string) – The filename

  • format (string) – Name of the format

Return type:

SequenceHandle

For files with non-standard extensions, the format can be set explicitly specifying the format parameter.

# recognizes FASTA file by file extension
myseq = io.LoadSequence('seq.fasta')
# for obtaining a SequenceList
seqlist = io.LoadSequenceList('seqs.fasta')
# or for multiple alignments (here from CLUSTAL)
aln = io.LoadAlignment('algnm.aln', format="clustal")

For a list of file formats supported by LoadSequence() see Supported Sequence File Formats.

Raises:

IOUnknownFormatException if the format string supplied is not recognized or the file format can not be detected based on the file extension.

IOException if the import fails due to an erroneous or inexistent file.

LoadSequenceList(filename, format='auto')

For a description of how to use LoadSequenceList() please refer to LoadSequence(). For a list of file formats supported by LoadSequenceList() see Supported Sequence File Formats.

Parameters:
  • filename (string) – The filename

  • format (string) – Name of the format

Return type:

SequenceList

LoadAlignment(filename, format='auto')

For a description of how to use LoadAlignment() please refer to LoadSequence(). For a list of file formats supported by LoadAlignment() see Supported Sequence File Formats.

Parameters:
  • filename (string) – The filename

  • format (string) – Name of the format

Return type:

AlignmentHandle

LoadSequenceProfile(filename, format='auto')

Load sequence profile data from disk. If format is set to ‘auto’, the function guesses the filetype based on the extension of the file. Files ending in ‘.hhm’ (output of HHblits) and ‘.pssm’ (ASCII Table (PSSM) output of PSI-BLAST as generated with blastpgp and flag -Q) will automatically be loaded.

For files with non-standard extensions, the format can be set explicitly specifying the format parameter.

# recognizes hhm file by file extension
myprof = io.LoadSequenceProfile('myhmm.hhm')
# recognizes pssm file by file extension
myprof = io.LoadSequenceProfile('myprof.pssm')

# to override format
myprof = io.LoadSequenceProfile('myfile', format='hhm')
myprof = io.LoadSequenceProfile('myfile', format='pssm')

For a list of file formats supported by LoadSequenceProfile() see Supported Sequence Profile File Formats.

Parameters:
  • filename (string) – The filename

  • format (string) – Name of the format

Return type:

ProfileHandle

Raises:

IOUnknownFormatException if the format string supplied is not recognized or the file format can not be detected based on the file extension.

IOException if the import fails due to an erroneous or inexistent file.

AlignmentFromString(data, format)

Load alignment from string.

The format argument is mandatory. For a list of supported formats, see Supported Sequence File Formats.

Parameters:
  • data (string) – The alignment

  • format (string) – Name of the format

Return type:

AlignmentHandle

SequenceFromString(data, format)

Load sequence from string.

The format argument is mandatory. For a list of supported formats, see Supported Sequence File Formats.

Parameters:
  • data (string) – The sequence

  • format (string) – Name of the format

Return type:

SequenceHandle

SequenceListFromString(data, format)

Load a list of sequences from string.

The format argument is mandatory. For a list of supported formats, see Supported Sequence File Formats.

Parameters:
  • data (string) – The list of sequences

  • format (string) – Name of the format

Return type:

SequenceList

SequenceProfileFromString(data, format)

Load sequence profile from string.

The format argument is mandatory.

Parameters:
  • data (string) – String containing the data you would read from a file with specified format.

  • format (string) – Name of the format. Can either be “hhm” or “pssm”.

Return type:

ProfileHandle

Saving Sequence Data

SaveSequence(sequence, filename, format='auto')

Saving sequence data is performed by calling SaveSequence(). For files with non-standard extensions, the format can be set explicitly specifying the ‘format’ parameter.

# recognizes FASTA file by file extension
io.SaveSequence(myseq, 'seq.fasta')
# for saving a SequenceList
io.SaveSequenceList(seqlist, 'seqlist.fasta')
# or for multiple alignments (here in FASTA format)
io.SaveAlignment(aln, 'aln.fasta')

For a list of file formats supported by SaveSequence() see Supported Sequence File Formats.

Parameters:
  • sequence (SequenceHandle) – The sequence

  • filename (string) – The filename

  • format (string) – Name of the format

Raises:

IOUnknownFormatException if the format string supplied is not recognized or the file format can not be detected based on the file extension.

IOException if the import fails due to an erroneous or inexistent file.

SaveSequenceList(seq_list, filename, format='auto')

For a desription of how to use SaveSequenceList() please refer to SaveSequence(). For a list of file formats supported by SaveSequenceList() see Supported Sequence File Formats.

Parameters:
  • seq_list (SequenceList) – The sequence list

  • filename (string) – The filename

  • format (string) – Name of the format

SaveAlignment(aln, filename, format='auto')

For a desription of how to use SaveAlignment() please refer to SaveSequence().

For a list of file formats supported by SaveAlignment() see Supported Sequence File Formats.

Parameters:
  • aln (AlignmentHandle) – The alignment

  • filename (string) – The filename

  • format (string) – Name of the format

AlignmentToString(ali, format)

Return alignment as a string.

The format argument is mandatory. For a list of supported formats see Supported Sequence File Formats.

Parameters:
  • ali (AlignmentHandle) – The alignment

  • format (string) – Name of the format

Return type:

string

SequenceToString(seq, format)

Return sequence as a string.

The format argument is mandatory. For a list of supported formats see Supported Sequence File Formats.

Parameters:
  • seq (SequenceHandle) – The sequence

  • format (string) – Name of the format

Return type:

string

SequenceListToString(seq_list, format)

Return sequence list as a string.

The format argument is mandatory. For a list of supported formats see Supported Sequence File Formats.

Parameters:
  • seq_list – The sequence list

  • format (string) – Name of the format

Return type:

string

Density Maps

Loading Density Maps

LoadImage(filename)
LoadImage(filename, format)

Load density map from disk. If no format is given, the function guesses the filetype based on the extension of the file. If the extension is unknown or not present the filetype will be guessed based on the content of the file if possible.

Parameters:
  • filename (string) – The filename

  • format – The file format

Raises:

IOUnknownFormatException if the format supplied is not recognized or the file format can not be detected based on the file extension and content.

IOException if the import fails due to an erroneous or inexistent file.

# recognizes mrc file by file extension
ent = io.LoadImage('file.mrc')

# it is always possible to explicitly set the image format
# DAT file explicitly
ent = io.LoadImage('file', Dat())

For a list of file formats supported by LoadImage(), see Supported Image File Formats.

Saving Density Maps

SaveImage(image, filename)
SaveImage(image, filename, format)

Save density map to disk. If no format is set, the function guesses the filetype based on the file extension.

Parameters:
  • image (IMageHandle) – The density map to be saved

  • filename (string) – The filename

  • format – The file format

Raises:

IOUnknownFormatException if the file format can not be detected based on the file extension

For a list of file formats supported by SaveImage(), see Supported Image File Formats.

# load density map
image = io.LoadImage('density_map.ccp4')
# save density map
io.SaveImage(image, 'new_map.map', CCP4())

Stereochemical Parameters

In order to check the structure for some stereo-chemical and steric clashes before computing the lDDT scores it is required to pass parameter file based on Engh and Huber parameters, and on the atomic radii as defined in the Cambridge Structural Database. OpenStructure ships with default file called stereo_chemical_props.txt located in $OST_ROOT/share/openstructure directory. A function ReadStereoChemicalPropsFile() is used to read this file.

ReadStereoChemicalPropsFile(filename='', check=True)

Read stereochemical parameters - if not provided a local version will be used.

Parameters:
  • filename (str) – The path to the parameter file that will be used. If set to “”, it reads the default file shipped with OpenStructure.

  • check (bool) – Raise an error when any of the resulting tables are empty.

Returns:

Object containing stereochemical parameters

Return type:

StereoChemicalProps

GetStereoChemicalPropsFile()

Get the default path to the stereochemical paramteres file.

Search

Enter search terms or a module, class or function name.

Contents

Documentation is available for the following OpenStructure versions:

(Currently viewing dev) / 2.7 / 2.6 / 2.5 / 2.4 / 2.3.1 / 2.3 / 2.2 / 2.1 / 2.0 / 1.9 / 1.8 / 1.7.1 / 1.7 / 1.6 / 1.5 / 1.4 / 1.3 / 1.2 / 1.11 / 1.10 / 1.1

This documentation is still under heavy development!
If something is missing or if you need the C++ API description in doxygen style, check our old documentation for further information.