chem

package module
v0.7.1 Latest Latest
Warning

This package is not in the latest version of its module.

Go to latest
Published: Apr 3, 2024 License: LGPL-2.1 Imports: 13 Imported by: 23

README

goChem is a library fro Computational Chemistry and Biochemistry written in the Go programming language.

Check out www.gochem.org for more information.

Design goals:

  • Simplicity.
  • Readability.
  • Fast and light
  • Concurrent when possible and necessary
  • Easy to extend
  • Useful for Computational Chemistry/Biochemistry at a classical and QM levels.

See the Wiki por goChem's current capabilities (https://github.com/rmera/gochem/wiki)

There is no publication on goChem yet. If you use goChem, or or any program based on goChem (such as goMD and goAnalize) before such a publication is available, we ask you to support goChem by citing the library in your publication as:

Domínguez, C., Jiménez, V, Savasci, G., Araya-Osorio, R., Pesonen, J., Mera-Adasme,  "goChem: A library for computational chemistry". https://www.gochem.org

Currently, gochem is licensed under LGPL2.1. This might change towards BSD in the future. Meanwhile, if you want to use some of this for a BSD-licensed project, contact the developer.

goChem uses gonum (https://www.gonum.org) for matrix operations. The user can choose to have gonum backed by a pure-go (and some assembly) implementation of BLAS (which allows not having runtime dependencies) or by the somewhat more efficient cBLAS.

Reading xtc files requires the xdrfile library from Gromacs (www.gromacs.org)

All dependencies of goChem are open source.

LICENSE

Copyright 2012 Raul Mera rmeraaatacademicosdotutadotcl

This program, including its documentation, is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2.1 of the License, or (at your option) any later version.

This program and its documentation is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with this program. If not, see http://www.gnu.org/licenses/.

The mascot is a modification by Sebastian Franchini of the Go language mascot by Renee French.

To the long life of the Ven. Khenpo Phuntzok Tenzin Rinpoche.

Documentation

Overview

This is the main package of the goChem library. It provides atom and molecule structures, facilities for reading and writing some files used in computational chemistry and functions for geometric manipulations and shape, among others indicators.

See www.gochem.org for more information.

  • ramacalc.go, part of gochem. * *
  • Copyright 2012 Raul Mera <rmera{at}chemDOThelsinkiDOTfi> *
  • This program is free software; you can redistribute it and/or modify
  • it under the terms of the GNU Lesser General Public License as
  • published by the Free Software Foundation; either version 2.1 of the
  • License, or (at your option) any later version. *
  • This program is distributed in the hope that it will be useful,
  • but WITHOUT ANY WARRANTY; without even the implied warranty of
  • MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  • GNU General Public License for more details. *
  • You should have received a copy of the GNU Lesser General
  • Public License along with this program. If not, see
  • <http://www.gnu.org/licenses/>. * *
  • Gochem is developed at the laboratory for instruction in Swedish, Department of Chemistry,
  • University of Helsinki, Finland. * *

Index

Constants

View Source
const (
	ErrNilData          = PanicMsg("goChem: Nil data given ")
	ErrInconsistentData = PanicMsg("goChem: Inconsistent data length ")
	ErrNilMatrix        = PanicMsg("goChem: Attempted to access nil v3.Matrix or gonum/mat64.Dense")
	ErrNilAtoms         = PanicMsg("goChem: Topology has a nil []*Atom slice")
	ErrNilAtom          = PanicMsg("goChem: Attempted to copy from or to a nil Atom")
	ErrAtomOutOfRange   = PanicMsg("goChem: Requested/Attempted setting Atom out of range")
	ErrNilFrame         = PanicMsg("goChem: Attempted to acces nil frame")
	ErrNotXx3Matrix     = PanicMsg("goChem: A v3.Matrix should have 3 columns")
	ErrCliffordRotation = PanicMsg("goChem-Clifford: Target and Result matrices must have the same dimensions. They cannot reference the same matrix") //the only panic that the Clifford functions throw.
)
View Source
const (
	Deg2Rad = 0.0174533
	Rad2Deg = 1 / Deg2Rad
	H2Kcal  = 627.509 //HArtree 2 Kcal/mol
	Kcal2H  = 1 / H2Kcal
	Kcal2KJ = 4.184
	KJ2Kcal = 1 / Kcal2KJ
	A2Bohr  = 1.889725989
	Bohr2A  = 1 / A2Bohr
	EV2Kcal = 23.061
	Kcal2EV = 1 / EV2Kcal
)

Conversions

View Source
const (
	CHDist = 1.098        //C(sp3)--H distance in A
	KBkJ   = 1.380649e-26 // Boltzmann constant kJ/K
	KB     = KBkJ * KJ2Kcal
	RkJ    = 8.31446261815324e-3 // kJ/(K*mol)
	NA     = 6.02214076e+23
	R      = RkJ * KJ2Kcal //.9872042586408e-3 //kcal/mol
)

Others

Variables

This section is empty.

Functions

func Angle

func Angle(v1, v2 *v3.Matrix) float64

Angle takes 2 vectors and calculate the angle in radians between them It does not check for correctness or return errors!

func AntiProjection

func AntiProjection(test, ref *v3.Matrix) *v3.Matrix

AntiProjection returns a vector in the direction of ref with the magnitude of a vector A would have if |test| was the magnitude of its projection in the direction of test.

func BestPlane

func BestPlane(coords *v3.Matrix, mol ...Masser) (*v3.Matrix, error)

BestPlane returns a row vector that is normal to the plane that best contains the molecule if passed a nil Masser, it will simply set all masses to 1. If more than one Masser is passed Only the first will be considered

func BestPlaneP

func BestPlaneP(evecs *v3.Matrix) (*v3.Matrix, error)

BestPlaneP takes sorted evecs, according to the eval,s and returns a row vector that is normal to the Plane that best contains the molecule. Notice that the function can't possibly check that the vectors are sorted. The P at the end of the name is for Performance. If That is not an issue it is safer to use the BestPlane function that wraps this one.

func BondedPaths added in v0.6.3

func BondedPaths(at *Atom, targetIndex int, options ...*BondedOptions) [][]int

BondedPaths determines the paths between at and the atom with Index targetIndex. It returns a slice of slices of int, where each sub-slice contains all the atoms in between at and the target (including the index of at) If there is no valid path, it returns nil. In the options structure BondeOption, OnlyShortest set to true causes only the shortest of the found paths to be returned. All atoms in the molecule need to have the "index" field filled. If the targetIndex is the same as that of the current atom, and the path is not given, nil, or of len 0, the function will search for a cyclic path back to the initial atom. if onlyshortest is true, only the shortest path will be returned (the other elements of the slice will be nil) This can be useful if you want to save memory on a very intrincate molecule.

func BondedPathsFunc added in v0.7.0

func BondedPathsFunc(at *Atom, targetIndex int, f func(*Bond) bool, options ...*BondedOptions) [][]int

If this works, we could replace BondedPaths by just a call to this function with f=func(*Bond)bool{return true} BondedPaths determines the paths between at and the atom with Index targetIndex. It returns a slice of slices of int, where each sub-slice contains all the atoms in between at and the target (including the index of at) If there is no valid path, it returns nil. In the options structure BondeOption, OnlyShortest set to true causes only the shortest of the found paths to be returned. All atoms in the molecule need to have the "index" field filled. If the targetIndex is the same as that of the current atom, and the path is not given, nil, or of len 0, the function will search for a cyclic path back to the initial atom. if onlyshortest is true, only the shortest path will be returned (the other elements of the slice will be nil) This can be useful if you want to save memory on a very intrincate molecule.

func CenterOfMass

func CenterOfMass(geometry *v3.Matrix, massS ...*mat.Dense) (*v3.Matrix, error)

CenterOfMass returns the center of mass the atoms represented by the coordinates in geometry and the masses in mass, and an error. If no mass is given, it calculates the geometric center

func Corrupted

func Corrupted(X Traj, R Atomer) error

Corrupted is a convenience function to check that a reference and a trajectory have the same number of atoms

func CutAlphaRef

func CutAlphaRef(r Atomer, chain []string, list []int) []int

CutAlphaRef will return a list with the atoms in the residues indicated by list, in the chains given. The carbonyl carbon and amide nitrogen for each residue will be transformer into hydrogens. The MolID of the other backbone atoms will be set to -1 so they are no longer considered.

func CutBackRef

func CutBackRef(r Atomer, chains []string, list [][]int) ([]int, error)

CutBackRef takes a list of lists of residues and selects from r all atoms in each the list list[i] and belonging to the chain chain[i]. It caps the N and C terminal of each list with -COH for the N terminal and NH2 for C terminal. the residues on each sublist should be contiguous to each other. for instance, {6,7,8} is a valid sublist, {6,8,9} is not. This is NOT currently checked by the function!. It returns the list of kept atoms

func CutBetaRef

func CutBetaRef(r Atomer, chain []string, list []int) []int

CutLateralRef will return a list with the atom indexes of the lateral chains of the residues in list for each of these residues it will change the alpha carbon to oxygen and change the residue number of the rest of the backbone to -1.

func Dihedral

func Dihedral(a, b, c, d *v3.Matrix) float64

Dihedral calculates the dihedral between the points a, b, c, d, where the first plane is defined by abc and the second by bcd. The ouputs are defined between 0 and pi, so it shouldn't be used for Ramachandran angles. Use RamaDihedral for that.

func DihedralAlt added in v0.6.0

func DihedralAlt(a, b, c, d *v3.Matrix) float64

DihedralAlt calculates the dihedral between the points a, b, c, d, where the first plane is defined by abc and the second by bcd. It is exactly the same as Dihedral, only kept for API stability.

func DihedralRama added in v0.7.0

func DihedralRama(a, b, c, d *v3.Matrix) float64

Similar to the Dihedral funcion, obatins the dihedral angle between 2 points. The difference is that the outputs for this functions are defined between -pi and pi, while Dihedral's outputs are defined between 0 and pi.

func EasyShape

func EasyShape(coords *v3.Matrix, epsilon float64, mol ...Masser) (float64, float64, error)

EasyShape takes a matrix of coordinates, a value for epsilon (a number close to zero, the closer, the more strict the orthogonality requriements are) and an (optative) masser and returns two shape indicators based on the elipsoid of inertia (or it massless equivalent) a linear and circular distortion indicators, and an error or nil (in that order). if you give a negative number as epsilon, the default (quite strict) will be used.

func EulerRotateAbout

func EulerRotateAbout(coordsorig, ax1, ax2 *v3.Matrix, angle float64) (*v3.Matrix, error)

EulerRotateAbout uses Euler angles to rotate the coordinates in coordsorig around by angle radians around the axis given by the vector axis. It returns the rotated coordsorig, since the original is not affected. It seems more clunky than the RotateAbout, which uses Clifford algebra. I leave it for benchmark, mostly, and might remove it later. There is no test for this function!

func FixGromacsPDB

func FixGromacsPDB(mol Atomer)

FixGromacsPDB fixes the problem that Gromacs PDBs have when there are more than 10000 residues Gromacs simply restarts the numbering. Since solvents (where this is likely to happen) don't have chain ID in Gromacs, it's impossible to distinguish between the water 1 and the water 10001. FixGromacsPDB Adds a chain ID to the newly restrated residue that is the letter/symbol coming after the last seen chain ID in the constant allchains defined in this file. The current implementation does nothing if a chain ID is already defined, even if it is wrong (if 9999 and the following 0 residue have the same chain).

func FixNumbering

func FixNumbering(r Atomer)

FixNumbering will put the internal numbering+1 in the atoms and residue fields, so they match the current residues/atoms in the molecule

func GroFileWrite added in v0.6.0

func GroFileWrite(outname string, Coords []*v3.Matrix, mol Atomer) error

GoFileWrite writes the molecule described by mol and Coords into a file in the Gromacs gro format. If Coords has more than one elements, it will write a multi-state file.

func GroSnapWrite added in v0.6.0

func GroSnapWrite(coords *v3.Matrix, mol Atomer, out io.Writer) error

GroSnapWrite writes a single snapshot of a molecule to an io.Writer, in the Gro format.

func Improper added in v0.6.0

func Improper(a, b, c, d *v3.Matrix) float64

Improper calculates the improper dihedral between the points a, b,c,d as the angle between the plane defined by a,b,c and that defined by the plane bcd

func ImproperAlt added in v0.6.0

func ImproperAlt(a, b, c, d *v3.Matrix) float64

ImproperAlt calculates the improper dihedral between the points a, b,c,d as the angle between the plane defined by a,b,c and the cd vector

func InWhichRing added in v0.6.3

func InWhichRing(at *Atom, rings []*Ring) int

InWhichRing returns the index of the first ring found to which the at atom belongs, or -1 if the atom is not part of any ring.

func MakeWater

func MakeWater(a1, a2 *v3.Matrix, distance, angle float64, oxygen bool) *v3.Matrix

MakeWater Creates a water molecule at distance Angstroms from a2, in a direction that is angle radians from the axis defined by a1 and a2. Notice that the exact position of the water is not well defined when angle is not zero. One can always use the RotateAbout function to move the molecule to the desired location. If oxygen is true, the oxygen will be pointing to a2. Otherwise, one of the hydrogens will.

func MassCenter

func MassCenter(in, oref *v3.Matrix, massS ...*mat.Dense) (*v3.Matrix, *v3.Matrix, error)

MassCenter centers in in the center of mass of oref. Mass must be A column vector. Returns the centered matrix and the displacement matrix.

func MassCenterMem added in v0.7.0

func MassCenterMem(in, oref, ret *v3.Matrix, massS ...*mat.Dense) (*v3.Matrix, error)

MassCenterMem centers in in the center of mass of oref, putting the result in ret. Mass must be A column vector. Returns the centered matrix and the displacement matrix.

func MemPerAtomRMSD added in v0.7.0

func MemPerAtomRMSD(test, templa, ctest, ctempla, tmp *v3.Matrix, indexes ...[]int) ([]float64, error)

MemMSD calculates the MSDs between test and templa, considering only the atoms present in the slices of int slices indexes. If given. If only one set of indexes is given, it will be assumed to beling to test, if it has less elements than that system, or to templa,otherwise.

func MemRMSD

func MemRMSD(test, templa, tmp *v3.Matrix, indexes ...[]int) (float64, error)

MemRMSD calculates the RMSD between test and template, considering only the atoms present in the slices of int slices indexes. The first indexes slices will be assumed to contain test indexes and the second, template indexes. If you give only one (it must be the first one), it will be assumed to correspond to whatever molecule that has more atoms than the elements in the slice. Giving a nil or 0-lenght first slice and a non-nil second slice will cause MemRMSD to not consider neither of them. The same number of atoms has to be considered for the calculation in both systems. It does not superimpose the objects. To save memory, it asks for the temporary matrix it needs to be supplied: tmp must be Nx3 where N is the number of elements in testlst and templalst NOTE: This function should ask for 3 tmp matrices, not just 1

func MolIDNameChain2Index

func MolIDNameChain2Index(mol Atomer, molID int, name, chain string) (int, error)

MolIDNameChain2Index takes a molID (residue number), atom name, chain index and a molecule Atomer. it returns the index associated with the atom in question in the Ref. The function returns also an error (if failure of warning) or nil (if succses and no warnings). Note that this function is not efficient to call several times to retrieve many atoms.

func Molecules2Atoms

func Molecules2Atoms(mol Atomer, residues []int, chains []string) []int

Molecules2Atoms gets a selection list from a list of residues. It select all the atoms that form part of the residues in the list. It doesnt return errors. If a residue is out of range, no atom will be returned for it. Atoms are also required to be part of one of the chains specified in chains, but a nil "chains" can be given to select all chains.

func MomentTensor

func MomentTensor(A *v3.Matrix, massslice ...[]float64) (*v3.Matrix, error)

MomentTensor returns the moment tensor for a matrix A of coordinates and a column vector mass with the respective massess.

func MultiPDBWrite

func MultiPDBWrite(out io.Writer, Coords []*v3.Matrix, mol Atomer, Bfactors [][]float64) error

MultiPDBWrite writes a multiPDB for the molecule mol and the various coordinate sets in CandB, to the io.Writer given. CandB is a list of lists of *matrix.DenseMatrix. If it has 2 elements or more, the second will be used as Bfactors. If it has one element, all b-factors will be zero. Returns an error if fails, or nil if succeeds.

func NegateIndexes added in v0.6.0

func NegateIndexes(indexes []int, length int) []int

NegateIndexes, given a set of indexes and the length of a molecule, produces a set of all the indexes _not_ in the original set.

func OnesMass

func OnesMass(lenght int) *v3.Matrix

OnesMass returns a column matrix with lenght rosw. This matrix can be used as a dummy mass matrix for geometric calculations.

func PDBFileWrite

func PDBFileWrite(pdbname string, coords *v3.Matrix, mol Atomer, Bfactors []float64) error

PDBFileWrite writes a PDB for the molecule mol and the coordinates Coords to a file name pdbname.

func PDBStringWrite

func PDBStringWrite(coords *v3.Matrix, mol Atomer, bfact []float64) (string, error)

PDBStringWrite writes a string in PDB format for a given reference, coordinate set and bfactor set, which must match each other returns the written string and error or nil.

func PDBWrite

func PDBWrite(out io.Writer, coords *v3.Matrix, mol Atomer, bfact []float64) error

PDBWrite writes a PDB formatted sequence of bytes to an io.Writer for a given reference, coordinate set and bfactor set, which must match each other. Returns error or nil.

func Projection

func Projection(test, ref *v3.Matrix) *v3.Matrix

Projection returns the projection of test in ref.

func RMSD

func RMSD(test, templa *v3.Matrix, indexes ...[]int) (float64, error)

RMSD calculates the RMSD between test and template, considering only the atoms present in the slices of int slices indexes. The first indexes slices will be assumed to contain test indexes and the second, template indexes. If you give only one, it will be assumed to correspondo to whatever molecule that has more atoms than the elements in the slice. The same number of atoms has to be considered for superposition in both systems. The objects are not superimposed before the calculation.

func RamaCalc

func RamaCalc(M *v3.Matrix, dihedrals []RamaSet) ([][]float64, error)

RamaCalc Obtains the values for the phi and psi dihedrals indicated in []Ramaset, for the structure M. The angles are in *degrees*. It returns a slice of 2-element slices, one for the phi the next for the psi dihedral, a and an error or nil.

func RhoShapeIndexes

func RhoShapeIndexes(rhos []float64) (float64, float64, error)

RhoShapeIndexes Get shape indices based on the axes of the elipsoid of inertia. linear and circular distortion, in that order, and error or nil. Based on the work of Taylor et al., .(1983), J Mol Graph, 1, 30 This function has NOT been tested thoroughly in the sense of the appropiateness of the indexes definitions.

func Rhos

func Rhos(momentTensor *v3.Matrix, epsilon ...float64) ([]float64, error)

Rhos returns the semiaxis of the elipoid of inertia given the the moment of inertia tensor.

func Rotate

func Rotate(Target, Res, axis *v3.Matrix, angle float64) *v3.Matrix

Rotate takes the matrix Target and uses Clifford algebra to _concurrently_ rotate each of its rows by angle radians around axis. Axis must be a 3D row vector. Target must be an N,3 matrix. The result is put in Res, which is also returned. This is a low-level function. Most commonly, you'll want to use RotateAbout instead.

func RotateAbout

func RotateAbout(coordsorig, ax1, ax2 *v3.Matrix, angle float64) (*v3.Matrix, error)

RotateAbout about rotates the coordinates in coordsorig around by angle radians around the axis given by the vector axis. It returns the rotated coordsorig, since the original is not affected. Uses Clifford algebra.

func RotateP

func RotateP(Target, Res, axis, Rv, Rvrev, tmp1, tmp2, tmp3, tmp4, itmp1, itmp2, itmp3 *v3.Matrix, angle float64, gorut int)

RotateP takes the matrix Target and uses Clifford algebra to _concurrently_ rotate each of its rows by angle radians around axis. Axis must be a 3D row vector. Target must be an N,3 matrix.

func RotateSer

func RotateSer(Target, ToRot, axis *v3.Matrix, angle float64) *v3.Matrix

RotateSer takes the matrix Target and uses Clifford algebra to rotate each of its rows by angle radians around axis. Axis must be a 3D row vector. Target must be an N,3 matrix. The Ser in the name is from "serial". ToRot will be overwritten and returned

func RotateSerP

func RotateSerP(Target, ToRot, axis, tmpv, Rv, Rvrev, Rotatedv, itmp1, itmp2, itmp3, itmp4, itmp5, itmp6 *v3.Matrix, angle float64)

RotateSerP takes the matrix Target and uses Clifford algebra to rotate each of its rows by angle radians around axis. Axis must be a 3D row vector. Target must be an N,3 matrix. The Ser in the name is from "serial". ToRot will be overwritten and returned. RotateSerP only allocates some floats but not any v3.Matrix. Instead, it takes the needed intermediates as arguments, hence the "P" for "performance" If performance is not an issue, use RotateSer instead, it will perform the allocations for you and call this function. Notice that if you use this function directly you may have to zero at least some of the intermediates before reusing them.

func RotatorAroundZ

func RotatorAroundZ(gamma float64) (*v3.Matrix, error)

RotatorAroundZ returns an operator that will rotate a set of coordinates by gamma radians around the z axis.

func RotatorAroundZToNewY

func RotatorAroundZToNewY(newy *v3.Matrix) (*v3.Matrix, error)

RotatorAroundZToNewY takes a set of coordinates (mol) and a vector (y). It returns a rotation matrix that, when applied to mol, will rotate it around the Z axis in such a way that the projection of newy in the XY plane will be aligned with the Y axis.

func RotatorToNewZ

func RotatorToNewZ(newz *v3.Matrix) *v3.Matrix

RotatorToNewZ takes a matrix a row vector (newz). It returns a linear operator such that, when applied to a matrix mol ( with the operator on the right side) it will rotate mol such that the z axis is aligned with newz.

func RotatorTranslatorToSuper

func RotatorTranslatorToSuper(test, templa *v3.Matrix) (*v3.Matrix, *v3.Matrix, *v3.Matrix, *v3.Matrix, error)

RotatorTranslatorToSuper superimposes the set of cartesian coordinates given as the rows of the matrix test on the gnOnes of the rows of the matrix templa. Returns the transformed matrix, the rotation matrix, 2 translation row vectors For the superposition plus an error. In order to perform the superposition, without using the transformed the first translation vector has to be added first to the moving matrix, then the rotation must be performed and finally the second translation has to be added. This is a low level function, although one can use it directly since it returns the transformed matrix. The math for this function is by Prof. Veronica Jimenez-Curihual, UNAB, Chile.

func ScaleBond

func ScaleBond(C, H *v3.Matrix, bond float64)

ScaleBond takes a C-H bond and moves the H (in place) so the distance between them is the one given (bond). CAUTION: I have only tested it for the case where the original distance>bond, although I expect it to also work in the other case.

func ScaleBonds

func ScaleBonds(coords *v3.Matrix, mol Atomer, n1, n2 string, finallenght float64)

ScaleBonds scales all bonds between atoms in the same residue with names n1, n2 to a final lenght finallengt, by moving the atoms n2. the o<ration is executed in place.

func SelCone

func SelCone(B, selection *v3.Matrix, angle, distance, thickness, initial float64, whatcone int) []int

SelCone, Given a set of cartesian points in sellist, obtains a vector "plane" normal to the best plane passing through the points. It selects atoms from the set A that are inside a cone in the direction of "plane" that starts from the geometric center of the cartesian points, and has an angle of angle (radians), up to a distance distance. The cone is approximated by a set of radius-increasing cilinders with height thickness. If one starts from one given point, 2 cgnOnes, one in each direction, are possible. If whatcone is 0, both cgnOnes are considered. if whatcone<0, only the cone opposite to the plane vector direction. If whatcone>0, only the cone in the plane vector direction. the 'initial' argument allows the construction of a truncate cone with a radius of initial.

func Super

func Super(test, templa *v3.Matrix, indexes ...[]int) (*v3.Matrix, error)

Super determines the best rotation and translations to superimpose the coords in test considering only the atoms present in the slices of int slices indexes. The first indexes slices will be assumed to contain test indexes and the second, template indexes. If you give only one, it will be assumed to correspond to test, if test has more atoms than elements on the indexes set, or templa, otherwise. If no indexes are given, all atoms on each system will be superimposed. The number of atoms superimposed on both systems must be equal. Super modifies the test matrix, but template and indexes are not touched.

func TagAtomsByName

func TagAtomsByName(r Atomer, name string, list []int) int

TagAtomsByName will tag all atoms with a given name in a given list of atoms. return the number of tagged atoms

func XYZFileAsTraj added in v0.6.0

func XYZFileAsTraj(xyzname string) (*Molecule, *XYZTraj, error)

Reads a multi-xyz file. Returns the first snapshot as a molecule, and the other ones as a XYZTraj

func XYZFileWrite

func XYZFileWrite(xyzname string, Coords *v3.Matrix, mol Atomer) error

XYZWrite writes the mol Ref and the Coord coordinates in an XYZ file with name xyzname which will be created fot that. If the file exist it will be overwritten.

func XYZStringWrite

func XYZStringWrite(Coords *v3.Matrix, mol Atomer) (string, error)

XYZStringWrite writes the mol Ref and the Coord coordinates in an XYZ-formatted string.

func XYZWrite

func XYZWrite(out io.Writer, Coords *v3.Matrix, mol Atomer) error

XYZWrite writes the mol Ref and the Coords coordinates to a io.Writer, in the XYZ format.

Types

type Atom

type Atom struct {
	Name string //PDB name of the atom
	ID   int    //The PDB index of the atom

	Tag       int     //Just added this for something that someone might want to keep that is not a float.
	MolName   string  //PDB name of the residue or molecule (3-letter code for residues)
	MolName1  byte    //the one letter name for residues and nucleotids
	Char16    byte    //Whatever is in the column 16 (counting from 0) in a PDB file, anything.
	MolID     int     //PDB index of the corresponding residue or molecule
	Chain     string  //One-character PDB name for a chain.
	Mass      float64 //hopefully all these float64 are not too much memory
	Occupancy float64 //a PDB crystallographic field, often used to store values of interest.
	Vdw       float64 //radius
	Charge    float64 //Partial charge on an atom
	Symbol    string
	Het       bool    // is the atom an hetatm in the pdb file? (if applicable)
	Bonds     []*Bond //The bonds connecting the atom to others.
	// contains filtered or unexported fields
}

Atom contains the information to represent an atom, except for the coordinates, which will be in a separate *v3.Matrix and the b-factors, which are in a separate slice of float64.

func (*Atom) Copy

func (N *Atom) Copy(A *Atom)

Copy returns a copy of the Atom object. puts the copy into the

func (*Atom) Index added in v0.6.3

func (N *Atom) Index() int

Index returns the index of the atom

func (*Atom) SetIndex added in v0.7.0

func (N *Atom) SetIndex(i int)

Index returns the index of the atom

type AtomMultiCharger

type AtomMultiCharger interface {
	Atomer

	//Charge gets the total charge of the topology
	Charge() int

	//Multi returns the multiplicity of the topology
	Multi() int
}

AtomChargerMultier is atomer but also gives a charge and multiplicity

type Atomer

type Atomer interface {

	//Atom returns the Atom corresponding to the index i
	//of the Atom slice in the Topology. Should panic if
	//out of range.
	Atom(i int) *Atom

	Len() int
}

Atomer is the basic interface for a topology.

type Bond added in v0.6.3

type Bond struct {
	Index  int
	At1    *Atom
	At2    *Atom
	Dist   float64
	Energy float64 //One might prefer to use energy insteaf of order.
	//NOTE
	//I'm not sure if I leave just the order and let the user "Abuse" that field for energy
	//or anything else they want, or if I keep this extra field.
	Order float64 //Order 0 means undetermined
}

Bond represents a chemical, covalent bond.

func (*Bond) Cross added in v0.6.3

func (B *Bond) Cross(origin *Atom) *Atom

Cross returns the atom bonded to the origin atom bond in the receiver.

func (*Bond) Remove added in v0.6.3

func (b *Bond) Remove() error

Remove removes the receiver bond from the the Bond slices in the corresponding atoms

type BondedOptions added in v0.6.3

type BondedOptions struct {
	OnlyShortest bool //Only return the shortest path between the atoms
	// contains filtered or unexported fields
}

BondedOptions contains options for the BondePaths function

type CError

type CError struct {
	// contains filtered or unexported fields
}

CError (Concrete Error) is the concrete error type for the chem package, that implements chem.Error

func (CError) Decorate

func (err CError) Decorate(dec string) []string

Decorate will add the dec string to the decoration slice of strings of the error, and return the resulting slice.

func (CError) Error

func (err CError) Error() string

type ConcTraj

type ConcTraj interface {

	//Is the trajectory ready to be read?
	Readable() bool

	/*NextConc takes a slice of bools and reads as many frames as elements the list has
	form the trajectory. The frames are discarted if the corresponding elemetn of the slice
	is false. The function returns a slice of channels through each of each of which
	a *matrix.DenseMatrix will be transmited*/
	NextConc(frames []*v3.Matrix) ([]chan *v3.Matrix, error)

	//Returns the number of atoms per frame
	Len() int
}

ConcTraj is an interface for a trajectory that can be read concurrently.

type Error

type Error interface {
	Error() string
	Decorate(string) []string //This is the new thing for errors. It allows you to add information when you pass it up. Each call also returns the "decoration" slice of strins resulting from the current call. If passed an empty string, it should just return the current value, not add the empty string to the slice.

}

Error is the interface for errors that all packages in this library implement. The Decorate method allows to add and retrieve info from the error, without changing it's type or wrapping it around something else.

type LastFrameError

type LastFrameError interface {
	TrajError
	NormalLastFrameTermination() //does nothing, just to separate this interface from other TrajError's

}

LastFrameError has a useless function to distinguish the harmless errors (i.e. last frame) so they can be filtered in a typeswith that looks for this interface.

type Masser

type Masser interface {

	//Returns a column vector with the massess of all atoms
	Masses() ([]float64, error)
}

Masser can return a slice with the masses of each atom in the reference.

type Molecule

type Molecule struct {
	*Topology
	Coords      []*v3.Matrix
	Bfactors    [][]float64
	XYZFileData []string //This can be anything. The main rationale for including it is that XYZ files have a "comment"
	// contains filtered or unexported fields
}

Molecule contains all the info for a molecule in many states. The info that is expected to change between states, Coordinates and b-factors are stored separately from other atomic info.

func GroFileRead added in v0.6.0

func GroFileRead(groname string) (*Molecule, error)

GroFileRead reads a file in the Gromacs gro format, returning a molecule.

func NewMolecule

func NewMolecule(coords []*v3.Matrix, ats Atomer, bfactors [][]float64) (*Molecule, error)

NewMolecule makes a molecule with ats atoms, coords coordinates, bfactors b-factors charge charge and unpaired unpaired electrons, and returns it. It doesnt check for consitency across slices or correct charge or unpaired electrons.

func PDBFileRead

func PDBFileRead(pdbname string, read_additional ...bool) (*Molecule, error)

PDBFileRead reads a pdb file from an io.Reader. Returns a Molecule. If there is one frame in the PDB the coordinates array will be of lenght 1. It also returns an error which is not really well set up right now. read_additional is now deprecated. The reader will just read

func PDBRead

func PDBRead(pdb io.Reader, read_additional ...bool) (*Molecule, error)

PDBRRead reads a pdb file from an io.Reader. Returns a Molecule. If there is one frame in the PDB the coordinates array will be of lenght 1. It also returns an error which is not really well set up right now. read_additional is now "deprecated", it will be set to true regardless. I have made it into

func Reduce

func Reduce(mol Atomer, coords *v3.Matrix, build int, report *os.File, executable string) (*Molecule, error)

Reduce uses the Reduce program (Word, et. al. (1999) J. Mol. Biol. 285, 1735-1747. For more information see http://kinemage.biochem.duke.edu) To protonate a protein and flip residues. It writes the report from Reduce to a file called Reduce.err in the current dir. The Reduce executable can be given, in "executable", otherwise it will be assumed that it is a file called "reduce" and it is in the PATH.

func XYZFileRead

func XYZFileRead(xyzname string) (*Molecule, error)

XYZFileRead Reads an xyz or multixyz file (as produced by Turbomole). Returns a Molecule and error or nil.

func XYZRead

func XYZRead(xyzp io.Reader) (*Molecule, error)

XYZRead Reads an xyz or multixyz formatted bufio.Reader (as produced by Turbomole). Returns a Molecule and error or nil.

func (*Molecule) AddFrame

func (M *Molecule) AddFrame(newframe *v3.Matrix)

AddFrame akes a matrix of coordinates and appends them at the end of the Coords. It checks that the number of coordinates matches the number of atoms.

func (*Molecule) AddManyFrames

func (M *Molecule) AddManyFrames(newframes []*v3.Matrix)

AddManyFrames adds the array of matrices newfames to the molecule. It checks that the number of coordinates matches the number of atoms.

func (*Molecule) Close added in v0.7.0

func (M *Molecule) Close()

Close just sets the "current" counter to 0. If you are using it as a trajectory, you can always just discard the molecule and let the CG take care of it, as there is nothing on disk linked to it..

func (*Molecule) Coord

func (M *Molecule) Coord(atom, frame int) *v3.Matrix

Coord returns the coords for the atom atom in the frame frame. panics if frame or coords are out of range.

func (*Molecule) Copy

func (M *Molecule) Copy(A *Molecule)

Copy puts in the receiver a copy of the molecule A including coordinates

func (*Molecule) Corrupted

func (M *Molecule) Corrupted() error

Corrupted checks whether the molecule is corrupted, i.e. the coordinates don't match the number of atoms. It also checks That the coordinate matrices have 3 columns.

func (*Molecule) Current

func (M *Molecule) Current() int

Current returns the number of the next read frame

func (*Molecule) Del

func (M *Molecule) Del(i int) error

Del Deletes atom i and its coordinates from the molecule.

func (*Molecule) DelCoord

func (M *Molecule) DelCoord(i int) error

DelCoord deletes the coodinate i from every frame of the molecule.

func (*Molecule) InitRead

func (M *Molecule) InitRead() error

InitRead initializes molecule to be read as a traj (not tested!)

func (*Molecule) Less

func (M *Molecule) Less(i, j int) bool

Less returns true if the value in the Bfactors for atom i are less than that for atom j, and false otherwise.

func (*Molecule) NFrames

func (M *Molecule) NFrames() int

NFrames returns the number of frames in the molecule

func (*Molecule) Next

func (M *Molecule) Next(V *v3.Matrix, box ...[]float64) error

Next puts the next frame into V and returns an error or nil The box argument is never used.

func (*Molecule) NextConc

func (M *Molecule) NextConc(frames []*v3.Matrix) ([]chan *v3.Matrix, error)

NextConc takes a slice of bools and reads as many frames as elements the list has form the trajectory. The frames are discarted if the corresponding elemetn of the slice is false. The function returns a slice of channels through each of each of which a *matrix.DenseMatrix will be transmited

func (*Molecule) Readable

func (M *Molecule) Readable() bool

Checks that the molecule exists and has some existent Coordinates, in which case returns true. It returns false otherwise. The coordinates could still be empty

func (*Molecule) SetCurrent

func (M *Molecule) SetCurrent(i int)

SetCurrent sets the value of the frame nex to be read to i.

func (*Molecule) Swap

func (M *Molecule) Swap(i, j int)

Swap function, as demanded by sort.Interface. It swaps atoms, coordinates (all frames) and bfactors of the molecule.

type PanicMsg

type PanicMsg string

PanicMsg is the type used for all the panics raised in the chem package so they can be easily recovered if needed. goChem only raises panics for programming errors: Attempts to to out of a matrice's bounds, dimension mismatches in matrices, etc.

func (PanicMsg) Error

func (v PanicMsg) Error() string

Error returns a string with an error message

type RamaSet

type RamaSet struct {
	Cprev   int
	N       int
	Ca      int
	C       int
	Npost   int
	MolID   int
	MolName string
}

A structure with the data for obtaining one pair of Ramachandran angles.

func RamaList

func RamaList(M Atomer, chains string, resran []int) ([]RamaSet, error)

RamaList takes a molecule and returns a slice of RamaSet, which contains the indexes for each dihedral to be included in a Ramachandran plot. It gets the dihedral indices for all residues in the range resran. if resran has 2 elements defining the boundaries. Otherwise, returns dihedral lists for the residues included in resran. If resran has 2 elements and the last is -1, RamaList will get all the dihedral for residues from resran[0] to the end of the chain. It only obtain dihedral lists for residues belonging to a chain included in chains.

func RamaResidueFilter

func RamaResidueFilter(dihedrals []RamaSet, filterdata []string, shouldBePresent bool) ([]RamaSet, []int)

RamaResidueFilter filters the set of a slice of RamaSet (i.e. a set of Ramachandran angles to be calculated) by residue.(ex. only GLY, everything but GLY) The 3 letter code of the residues to be filtered in or out is in filterdata, whether they are filter in or out depends on shouldBePresent. It returns the filtered data and a slice containing the indexes in the new data of the residues in the old data, when they are included, or -1 when they are not included.

type Ring added in v0.6.3

type Ring struct {
	Atoms []int
	// contains filtered or unexported fields
}

Ring represents a molecular cycle.

func FindRings added in v0.6.3

func FindRings(coords *v3.Matrix, mol Atomer) []*Ring

Identifies and returns all rings in mol, by searching for cyclic paths.

func (*Ring) AddHs added in v0.6.3

func (R *Ring) AddHs(mol Atomer)

AddHs Adds to the ring the H atoms bonded to its members mol is the _entire_ molecule of which the receiver ring is part. It will panic if an atom of the ring is not found in mol, or is not bonded to anything

func (*Ring) IsIn added in v0.6.3

func (R *Ring) IsIn(index int) bool

IsIn returns true or false depending on whether the atom with the given index is part of the ring

func (*Ring) Planarity added in v0.6.3

func (R *Ring) Planarity(coord *v3.Matrix) float64

Planarity returns the planarity percentage of the receiver ring coords is the set of coordinates for the _entire_ molecule of which the ring is part. Planarity does not check that coords indeed corresponds to the correct molecule, so, doing so is the user's responsibility.

func (*Ring) Size added in v0.6.3

func (R *Ring) Size() int

Size returns the number of atoms in the ring

type Topology

type Topology struct {
	Atoms []*Atom
	// contains filtered or unexported fields
}

Topology contains information about a molecule which is not expected to change in time (i.e. everything except for coordinates and b-factors)

func BackboneCGize added in v0.6.0

func BackboneCGize(coord *v3.Matrix, mol Atomer, top bool) (*v3.Matrix, *Topology, error)

BackboneCGize takes a coord and a mol for a protein, and returns a new set of coordinates each of which cooresponds to the center of mass of the backbone of the corresponding residue in the original molecule. If top is true, it also returns a topology where each atom corrsponds to the center of mass, with the name "BB", and the correct residue name and ID. Otherwise it returns an empty topology. In other words, it takes a protein and returns a CG model for its backbone, in the currect conformation.

func MergeAtomers

func MergeAtomers(A, B Atomer) *Topology

Merges A and B in a single topology which is returned

func NewTopology

func NewTopology(charge, multi int, ats ...[]*Atom) *Topology

NewTopology returns topology with ats atoms, charge charge and multi multiplicity. It doesnt check for consitency across slices, correct charge or unpaired electrons.

func (*Topology) AppendAtom

func (T *Topology) AppendAtom(at *Atom)

AppendAtom appends an atom at the end of the reference

func (*Topology) AssignBonds added in v0.6.3

func (T *Topology) AssignBonds(coord *v3.Matrix) error

AssignBonds assigns bonds to a molecule based on a simple distance criterium, similar to that described in DOI:10.1186/1758-2946-3-33

func (*Topology) Atom

func (T *Topology) Atom(i int) *Atom

Atom returns the Atom corresponding to the index i of the Atom slice in the Topology. Panics if out of range.

func (*Topology) Charge

func (T *Topology) Charge() int

Charge returns the total charge of the topology

func (*Topology) CopyAtoms

func (T *Topology) CopyAtoms(A Atomer)

CopyAtoms copies the atoms form A into the receiver topology. This is a deep copy, so the receiver must have at least as many atoms as A.

func (*Topology) DelAtom

func (T *Topology) DelAtom(i int)

DelAtom Deletes atom i by reslicing. This means that the copy still uses as much memory as the original T.

func (*Topology) FillIndexes added in v0.6.3

func (T *Topology) FillIndexes()

FillsIndexes sets the Index value of each atom to that cooresponding to its place in the molecule.

func (*Topology) FillMasses added in v0.6.3

func (T *Topology) FillMasses()

FillMasses tries to get fill the masses for atom that don't have one by getting it from the symbol. Only a few common elements are supported

func (*Topology) FillVdw added in v0.6.3

func (T *Topology) FillVdw()

FillVdw tries to get fill the van der Waals radii for the atoms in the molecule from a symbol->radii map. Only a few common elements are supported

func (*Topology) Len

func (T *Topology) Len() int

Len returns the number of atoms in the topology.

func (*Topology) Masses

func (T *Topology) Masses() ([]float64, error)

Masses returns a slice of float64 with the masses of the atoms in the topology, or nil and an error if they have not been calculated

func (*Topology) Multi

func (T *Topology) Multi() int

Multi returns the multiplicity in the topology

func (*Topology) ResetIDs

func (T *Topology) ResetIDs()

ResetIDs sets the current order of atoms as ID and the order of molecules as MolID for all atoms

func (*Topology) SetAtom

func (T *Topology) SetAtom(i int, at *Atom)

SetAtom sets the (i+1)th Atom of the topology to aT. Panics if out of range

func (*Topology) SetCharge

func (T *Topology) SetCharge(i int)

SetCharge sets the total charge of the topology to i

func (*Topology) SetMulti

func (T *Topology) SetMulti(i int)

SetMulti sets the multiplicity in the topology to i

func (*Topology) SomeAtoms

func (R *Topology) SomeAtoms(T Atomer, atomlist []int)

SelectAtoms puts the subset of atoms in T that have indexes in atomlist into the receiver. Panics if problem.

func (*Topology) SomeAtomsSafe

func (R *Topology) SomeAtomsSafe(T Atomer, atomlist []int) error

SelectAtomsSafe puts the atoms of T with indexes in atomlist into the receiver. Returns error if problem.

type Traj

type Traj interface {

	//Is the trajectory ready to be read?
	Readable() bool

	//reads the next frame and returns it as DenseMatrix if keep==true, or discards it if false
	//it can also fill the (optional) box with the box vectors, it present in the frame.
	Next(output *v3.Matrix, box ...[]float64) error

	//Returns the number of atoms per frame
	Len() int
}

Traj is an interface for any trajectory object, including a Molecule Object

type TrajError

type TrajError interface {
	Error
	Critical() bool
	FileName() string
	Format() string
}

TrajError is the nterface for errors in trajectories

type XYZTraj added in v0.6.0

type XYZTraj struct {
	// contains filtered or unexported fields
}

XYZTraj is a trajectory-like representation of an XYZ File.

func (*XYZTraj) Len added in v0.6.0

func (X *XYZTraj) Len() int

func (*XYZTraj) Next added in v0.6.0

func (X *XYZTraj) Next(coords *v3.Matrix, box ...[]float64) error

Next reads the next snapshot of the trajectory into coords, or discards it, if coords is nil. It can take a box slice of floats, but won't do anything with it (only for compatibility with the Traj interface.

func (*XYZTraj) Readable added in v0.6.0

func (X *XYZTraj) Readable() bool

Readable returns true if the trajectory is fit to be read, false otherwise.

Directories

Path Synopsis
traj
dcd
stf
xtc

Jump to

Keyboard shortcuts

? : This menu
/ : Search site
f or F : Jump to
y or Y : Canonical URL