LOOS 4.1.0
The Lightweight Object Oriented Structural analysis library/toolkit
Loading...
Searching...
No Matches
Namespaces | Classes | Typedefs | Functions | Variables
loos Namespace Reference

Namespace for most things not already encapsulated within a class. More...

Namespaces

namespace  DensityTools
 Namespace for Density package.
 
namespace  internal
 Loos esoterica.
 
namespace  Math
 Namespace for math and math-related things in loos.
 
namespace  OptionsFramework
 Namespace for encapsulating options processing.
 

Classes

class  AbstractObserver
 
class  Amber
 Class for reading in AMBER parmtop/coord files... More...
 
class  AmberNetcdf
 Class for reading Amber Trajectories in NetCDF format. More...
 
class  AmberRst
 Class for reading amber restart files as a single-frame trajectory. More...
 
class  AmberTraj
 Class for reading amber coordinate trajectories. More...
 
struct  AndSelector
 Combines two selectors with a logical "and". More...
 
class  AscendingSort
 Policy class for sorting in ascending sequence. More...
 
class  Atom
 Basic Atom class for handling atom properties. More...
 
struct  AtomCoordsEquals
 Compares two atoms based on name, id, resid, resname, segid, and coords. More...
 
struct  AtomEquals
 Compares two atoms based solely on name, id, resid, resname, and segid. More...
 
class  AtomicGroup
 Class for handling groups of Atoms (pAtoms, actually) More...
 
struct  AtomNameSelector
 Predicate for selecting atoms based on explicit name matching. More...
 
struct  AtomSelector
 Virtual base-class for selecting atoms from a group. More...
 
class  BackboneSelector
 Predicate for selecting backbone. More...
 
class  BasicProgress
 
struct  BoundFmt
 Internal helper class to bind formatting state. More...
 
struct  CAlphaSelector
 Predicate for selecting CA atoms. More...
 
class  CCPDB
 Class for interpreting concatenated PDB files as a Trajectory. More...
 
class  CHARMM
 Class for reading a CHARMM coordinate file. More...
 
class  Coord
 Basic 3-D coordinates class. More...
 
class  DCD
 Class for reading DCD files. More...
 
class  DCDWriter
 A very lightweight class for writing simple DCDs. More...
 
class  DescendingSort
 Policy class for sorting in descending sequence. More...
 
class  EstimatingCounter
 A progress counter that can estimate how much time is left. More...
 
class  FileError
 Errors related to File I/O. More...
 
class  FileOpenError
 Error while opening a file. More...
 
class  FileReadError
 Errors that occur while reading a file. More...
 
class  FileReadErrorWithLine
 Errors that occur while reading a text file (where lines are tracked) More...
 
class  FileWriteError
 Errors while writing to files. More...
 
class  Fmt
 Output formatter class, adapted from Stroustrup's book. More...
 
class  FormFactor
 
class  FormFactorSet
 
class  Gromacs
 Implements a GROMACS model file (.gro) More...
 
class  HBondDetector
 Class for detecting hydrogen bonds. More...
 
struct  HeavyAtomSelector
 Select non-hydrogen atoms. More...
 
struct  HeavySolventSelector
 Select only heavy solvent atoms. More...
 
struct  HydrogenSelector
 Select hydrogen atoms. More...
 
struct  IndexSelector
 Predicate for selecting atoms based on index. More...
 
class  Kernel
 The Kernel (virtual machine) for compiling and executing user-defined atom selections. More...
 
class  KernelSelector
 Selection predicate that executes a compiled Kernel. More...
 
class  LineReader
 Class for reading line-by-line from a file while tracking line numbers and stripping comments. More...
 
class  LOOSError
 Generic LOOS exception. More...
 
class  Matrix44
 Specialized 4x4 Matrix class for handling coordinate transforms. More...
 
class  MatrixReadError
 Generic reading error class. More...
 
struct  MatrixReadImpl
 
struct  MatrixReadImpl< T, Math::Triangular, S >
 Special handling for triangular matrices. More...
 
struct  MatrixReadImpl< T, P, Math::SparseArray >
 Special handling for sparse matrices. More...
 
struct  MatrixWriteImpl
 
struct  MatrixWriteImpl< T, Math::Triangular, S, F >
 Write out a triangular matrix. More...
 
struct  MatrixWriteImpl< T, P, Math::SparseArray, F >
 Write out a sparse matrix. More...
 
class  MDTraj
 Class for reading a MDTraj HDF5 file. More...
 
class  MDTrajTraj
 Class for reading MDTraj HDF5 coordinate trajectories. More...
 
class  MMCIF
 Class to read pdbx/mmcif files. More...
 
class  MultiTrajectory
 Combine multiple trajectories (pTraj's) into one large virtual trajectory. More...
 
struct  NotSelector
 Negates a selection predicate. More...
 
class  NumericalError
 Exception caused by a blas/atlas error. More...
 
class  OptionsError
 Exception in options. More...
 
struct  OrSelector
 Combines two selectors with a logical "or". More...
 
class  ParseError
 Exception when parsing input data. More...
 
class  Parser
 Front-end to the Bison/Flex parser. More...
 
struct  ParserDriver
 
class  PDB
 PDB reading/writing class. More...
 
class  PDBTraj
 Class for intepreting separate PDB files as a single trajectory. More...
 
class  PercentProgress
 Provide feedback by percent-complete. More...
 
class  PercentProgressWithTime
 
class  PercentTrigger
 Trigger whenever at least frac percent more iterations have happened. More...
 
class  PeriodicBox
 Class for managing periodic box information. More...
 
class  PreciseMatrixFormatter
 Generic matrix element formatter allowing setting of width and precision. More...
 
class  ProgressCounter
 The progress counter front-end. More...
 
class  PSF
 Class for reading a subset of the PSF format. More...
 
class  Remarks
 Class for handling PDB Remarks. More...
 
struct  ResidRangeSelector
 Predicate for selecting atoms from a range of resid's. More...
 
class  RnaSuite
 Class for assigning backbone suites to an RNA. More...
 
class  ScientificMatrixFormatter
 Generic matrix element formatter using scientific notation... More...
 
struct  SegidSelector
 Predicate for selecting atoms based on the passed segid string. More...
 
class  SharedPeriodicBox
 This class manages a shared Periodicbox. More...
 
class  SimpleCounter
 Basic progress counter object, defining the interface... More...
 
struct  SolventSelector
 Predicate for selecting solvent based on common solvent SEGIDs. More...
 
class  StreamWrapper
 Simple wrapper class for caching stream pointers. More...
 
class  Timer
 Class for tracking time. More...
 
class  TimeSeries
 Time Series Class. More...
 
class  TinkerArc
 Class for handling Tinker ARC files (concatenation of .xyz files) More...
 
class  TinkerXYZ
 Class for reading a subset of the TinkerXYZ format. More...
 
class  Trajectory
 Base-class for polymorphic trajectories. More...
 
class  TrajectoryError
 Errors related to trajectory reading and writing. More...
 
class  TrajectoryWriter
 Base class for writing trajectories. More...
 
class  TriggerAlways
 This is a simple "trigger" for use as a default. More...
 
class  TriggerEvery
 Trigger every i-iterations. More...
 
class  TRR
 Class representing the GROMACS TRR trajectory files. More...
 
class  UniformWeight
 
class  UniqueStrings
 Class for uniquifying strings... More...
 
class  UnitCell
 This class encapsulates crystallographic unit cell data. More...
 
struct  unordered_pair_eq
 test equality of unordered pair. typename T must have == operator. More...
 
struct  unordered_pair_hash
 
class  UnsetProperty
 Exception when trying to use an unset Atom property. More...
 
class  UserTimer
 Policy class for tracking only user process time. More...
 
class  WallTimer
 Policy class for tracking wall-time. More...
 
class  Weights
 
class  WeightsFromFile
 
struct  XDRDataSizeError
 Exception indicating internal XDR error. More...
 
class  XForm
 Matrix class for handling coordinate transforms... More...
 
class  XTC
 Class representing GROMACS reduced precision, compressed trajectories. More...
 
class  XTCWriter
 Class for writing Gromacs XTC trajectories. More...
 
struct  ZSliceSelector
 Predicate for selecting atoms in a specific range of z values. More...
 

Typedefs

typedef boost::unordered_map< int, int > IMap
 
typedef boost::shared_ptr< AtomicGrouppAtomicGroup
 
typedef double greal
 
typedef long gint
 
typedef float dcd_real
 
typedef double dcd_double
 
typedef Coord< double > GCoord
 
typedef boost::shared_ptr< GCoordpGCoord
 
typedef boost::shared_ptr< AtompAtom
 
typedef boost::shared_ptr< TrajectorypTraj
 
typedef boost::shared_ptr< DCDpDCD
 
typedef boost::shared_ptr< AmberTrajpAmberTraj
 
typedef boost::shared_ptr< AmberNetcdfpAmberNetcdf
 
typedef boost::shared_ptr< CCPDBpCCPDB
 
typedef boost::shared_ptr< TinkerArcpTinkerArc
 
typedef boost::shared_ptr< PDBTrajpPDBTraj
 
typedef boost::shared_ptr< XTCpXTC
 
typedef boost::shared_ptr< TRRpTRR
 
typedef boost::shared_ptr< TrajectoryWriterpTrajectoryWriter
 
typedef boost::shared_ptr< PDBpPDB
 
typedef boost::shared_ptr< PSFpPSF
 
typedef boost::shared_ptr< AmberpAmber
 
typedef boost::shared_ptr< AmberRstpAmberRst
 
typedef boost::shared_ptr< TinkerXYZpTinkerXYZ
 
typedef boost::shared_ptr< GromacspGromacs
 
typedef boost::shared_ptr< CHARMMpCHARMM
 
typedef boost::shared_ptr< MMCIFpMMCIF
 
typedef Math::Matrix< float, Math::ColMajorRealMatrix
 
typedef Math::Matrix< double, Math::ColMajorDoubleMatrix
 
typedef TimeSeries< double > dTimeSeries
 
typedef TimeSeries< float > fTimeSeries
 
typedef boost::mt19937 base_generator_type
 
typedef Matrix44< greal > GMatrix
 

Functions

boost::tuple< std::vector< XForm >, greal, int > iterativeAlignment (alignment::vecMatrix &ensemble, greal threshold, int maxiter)
 
boost::tuple< std::vector< XForm >, greal, int > iterativeAlignment (std::vector< AtomicGroup > &ensemble, greal threshold, int maxiter)
 
boost::tuple< std::vector< XForm >, greal, int > iterativeAlignment (const AtomicGroup &model, pTraj &traj, const std::vector< uint > &frame_indices, greal threshold=1e-6, int maxiter=1000)
 Compute an iterative superposition by reading in frames from the Trajectory.
 
boost::tuple< std::vector< XForm >, greal, int > iterativeAlignment (const AtomicGroup &g, pTraj &traj, greal threshold, int maxiter)
 
bool isFileNetCDF (const std::string &fname)
 Returns true if the file is a NetCDF file.
 
std::ostream & operator<< (std::ostream &os, const loos::Atom &a)
 
AtomicGroup operator+ (const pAtom &lhs, const pAtom &rhs)
 
AtomicGroup operator+ (const pAtom &lhs, const AtomicGroup &rhs)
 
std::ostream & operator<< (std::ostream &os, const AtomicGroup &grp)
 
unsigned int deduceAtomicNumberFromMass (const double mass, const double tolerance=0.1)
 Deduce an atomic number from the mass.
 
std::string deduceElementNameFromMass (const double mass, const double tolerance)
 
double deduceMassFromAtomicNumber (const unsigned int number)
 Deduce a mass from the atomic number.
 
template<class T >
Coord< T > operator* (const Matrix44< T > &, const Coord< T > &)
 
double * gridify (DCD &, double *, double *, int[], const std::vector< int >, int, int, double)
 
AtomicGroup averageStructure (const std::vector< AtomicGroup > &ensemble)
 Compute the average structure of a set of AtomicGroup objects.
 
AtomicGroup averageStructure (const std::vector< AtomicGroup > &ensemble, const std::vector< XForm > &xforms)
 Compute the average structure of a set of AtomicGroup objects.
 
AtomicGroup averageStructure (const AtomicGroup &, const std::vector< XForm > &, pTraj &traj, const std::vector< uint > &indices)
 Compute the average structure from a trajectory reading only certain frames.
 
AtomicGroup averageStructure (const AtomicGroup &, const std::vector< XForm > &, pTraj &traj)
 Compute the average structure using all frames in a trajectory.
 
void applyTransforms (std::vector< AtomicGroup > &ensemble, std::vector< XForm > &xforms)
 
void readTrajectory (std::vector< AtomicGroup > &ensemble, const AtomicGroup &model, pTraj trajectory)
 
void readTrajectory (std::vector< AtomicGroup > &ensemble, const AtomicGroup &model, pTraj trajectory, std::vector< uint > &frames)
 
RealMatrix extractCoords (const std::vector< AtomicGroup > &ensemble)
 
RealMatrix extractCoords (const std::vector< AtomicGroup > &ensemble, const std::vector< XForm > &xforms)
 
void subtractAverage (RealMatrix &M)
 
boost::tuple< RealMatrix, RealMatrix, RealMatrixsvd (std::vector< AtomicGroup > &ensemble, const bool align=true)
 Compute the SVD of an ensemble with optional alignment (note RSVs returned are transposed)
 
void appendCoords (std::vector< std::vector< double > > &M, AtomicGroup &model, pTraj &traj, const std::vector< uint > &indices, const bool updates=false)
 
std::vector< std::vector< double > > readCoords (AtomicGroup &model, pTraj &traj, const std::vector< uint > &indices, const bool updates=false)
 
std::ostream & operator<< (std::ostream &os, const BoundFmt &bf)
 Create the output with the specified formatter.
 
std::ostream & operator<< (std::ostream &os, const Gromacs &g)
 Output the group as a GRO...
 
vector< uintparseIndexRange (const std::string &input, const uint maxsize)
 Use boost::Spirit to parse a new-style range list.
 
std::ostream & operator<< (std::ostream &os, const Kernel &k)
 
template<class T >
std::ostream & operator<< (std::ostream &os, const Timer< T > &t)
 
template<class T , class P , template< typename > class S>
Math::Matrix< T, P, S > readAsciiMatrix (std::istream &is)
 Read in a matrix from a stream returning a newly created matrix.
 
template<class T , class P , template< typename > class S>
void readAsciiMatrix (std::istream &is, Math::Matrix< T, P, S > &M)
 Read in a matrix from a stream storing it in the specified matrix.
 
template<class T , class P , template< typename > class S>
Math::Matrix< T, P, S > readAsciiMatrix (const std::string &fname)
 Read in a matrix from a file returning a newly created matrix.
 
template<class T , class P , template< typename > class S>
void readAsciiMatrix (const std::string &fname, Math::Matrix< T, P, S > &M)
 Read in a matrix from a file storing it in the specified matrix.
 
template<class T , class P , template< typename > class S, class F >
std::ostream & writeAsciiMatrix (std::ostream &os, const Math::Matrix< T, P, S > &M, const std::string &meta, const Math::Range &start, const Math::Range &end, const bool trans=false, F fmt=F())
 Write a submatrix to a stream.
 
template<class T , class P , template< typename > class S>
std::ostream & writeAsciiMatrix (std::ostream &os, const Math::Matrix< T, P, S > &M, const std::string &meta, const Math::Range &start, const Math::Range &end, const bool trans=false)
 Write a submatrix to a stream.
 
template<class T , class P , template< typename > class S, class F >
std::ostream & writeAsciiMatrix (std::ostream &os, const Math::Matrix< T, P, S > &M, const std::string &meta, const bool trans=false, F fmt=F())
 Write an entire matrix to a stream.
 
template<class T , class P , template< typename > class S>
std::ostream & writeAsciiMatrix (std::ostream &os, const Math::Matrix< T, P, S > &M, const std::string &meta, const bool trans=false)
 Write an entire matrix to a stream.
 
template<class T , class P , template< typename > class S, class F >
void writeAsciiMatrix (const std::string &fname, const Math::Matrix< T, P, S > &M, const std::string &meta, const Math::Range &start, const Math::Range &end, const bool trans=false, F fmt=F())
 Write a submatrix to a file.
 
template<class T , class P , template< typename > class S>
void writeAsciiMatrix (const std::string &fname, const Math::Matrix< T, P, S > &M, const std::string &meta, const Math::Range &start, const Math::Range &end, const bool trans=false)
 Write a submatrix to a file.
 
template<class T , class P , template< typename > class S, class F >
void writeAsciiMatrix (const std::string &fname, const Math::Matrix< T, P, S > &M, const std::string &meta, const bool trans=false, F fmt=F())
 Write an entire matrix to a file.
 
template<class T , class P , template< typename > class S>
void writeAsciiMatrix (const std::string &fname, const Math::Matrix< T, P, S > &M, const std::string &meta, const bool trans=false)
 Write an entire matrix to a file.
 
std::ostream & FormattedUnitCell (std::ostream &os, const UnitCell &u)
 
std::ostream & XTALLine (std::ostream &os, const GCoord &box)
 
std::ostream & FormatConectRecords (std::ostream &os, const PDB &p)
 
std::ostream & operator<< (std::ostream &os, const PDB &p)
 Output the group as a PDB...
 
std::ostream & operator<< (std::ostream &os, const Remarks &r)
 
std::string availableSystemFileTypes (const std::string &prefix)
 
pAtomicGroup createSystemPtr (const std::string &filename, const std::string &filetype)
 
pAtomicGroup createSystemPtr (const std::string &filename)
 
AtomicGroup createSystem (const std::string &filename)
 Factory function for reading in structure files.
 
AtomicGroup createSystem (const std::string &filename, const std::string &filetype)
 
std::string availableTrajectoryFileTypes (const std::string &prefix)
 
pTraj createTrajectory (const std::string &filename, const std::string &filetype, const AtomicGroup &g)
 
pTraj createTrajectory (const std::string &, const AtomicGroup &)
 Factory function for reading in a trajectory file.
 
std::string availableOutputTrajectoryFileTypes (const std::string &prefix)
 
pTrajectoryWriter createOutputTrajectory (const std::string &filename, const std::string &type, const bool append)
 
pTrajectoryWriter createOutputTrajectory (const std::string &filename, const bool append)
 
template<typename T , class SortPredicate >
std::vector< uint > sortedIndex (const T &A)
 
template<typename T >
std::vector< uint > sortedIndex (const T &A)
 Sort a container in ascending sequence.
 
std::string findBaseName (const std::string &)
 Pull off the file name extension (if present)
 
boost::tuple< std::string, std::string > splitFilename (const std::string &filename)
 
std::string getNextLine (std::istream &is, int *lineno)
 Get the next line of input, skipping blanks and stripping comments.
 
std::string invocationHeader (int, char *[])
 Create an invocation header.
 
std::vector< int > parseRangeList (const std::string &, const int endpoint=0)
 Parses a list of Octave-style range specifiers (for compatability)
 
AtomicGroup selectAtoms (const AtomicGroup &, const std::string)
 Applies a string-based selection to an atomic group...
 
std::string timeAsString (const double t, const uint precision=0)
 Convert t (seconds) into a string, converting to hours and minutes as necessary.
 
template<>
std::string parseStringAs< std::string > (const std::string &source, const uint pos, const uint nelem)
 
template<>
std::string fixedSizeFormat (const std::string &s, const uint n)
 
int parseStringAsHybrid36 (const std::string &source, const uint pos=0, const uint nelem=0)
 Convert a hybrid-36 encoded string into an int.
 
std::string hybrid36AsString (int value, uint fieldsize)
 Convert an int into a hybrid-36 encoded string.
 
std::string sanitizeString (const std::string &s)
 Removes internal newlines from string.
 
std::string stringsAsComments (const std::vector< std::string > &v)
 Converts a vector of strings into a standard log format.
 
std::string stringsAsString (const std::vector< std::string > &v)
 Converts a vector of strings into a single string with newlines.
 
template<>
std::string vectorAsStringWithCommas (const std::vector< std::string > &v)
 Specialization for strings that sanitizes the contained strings.
 
long availableMemory ()
 
template<typename T >
std::vector< T > readVector (LineReader &reader)
 Read a list of items using a LineReader object.
 
template<typename T >
std::vector< T > readVector (std::istream &is)
 Read a list of items from a stream with default behavior.
 
template<typename T >
std::vector< T > readVector (const std::string &fname)
 Read a list of items from a file with default behavior.
 
template<typename T >
std::vector< std::vector< T > > readTable (LineReader &reader)
 Read in a table of items using a LineReader object.
 
template<typename T >
std::vector< std::vector< T > > readTable (std::istream &is)
 Read in a table given a stream.
 
template<typename T >
std::vector< std::vector< T > > readTable (const std::string &fname)
 Read in a table given a filename.
 
template<typename T >
std::vector< T > parseRange (const std::string &text, const T endpoint=0)
 Parse an Octave/Matlab-style range.
 
template<typename T >
std::vector< T > parseRangeList (const std::string &text, const T endpoint=0)
 Parses a comma-separated list of Octave-style ranges.
 
template<typename T >
std::vector< T > parseRangeList (const std::vector< std::string > &ranges, const T endpoint=0)
 Parses a list of Octave-style ranges taken from a vector of strings.
 
template<typename T >
swab (const T &datum)
 Returns a byte-swapped copy of an arbitrary type.
 
template<typename T >
parseStringAs (const std::string &source, const uint pos=0, const uint nelem=0)
 Extracts a field from a string.
 
template<>
std::string parseStringAs< std::string > (const std::string &source, const uint pos, const uint nelem)
 
template<typename T >
std::string fixedSizeFormat (const T t, const uint n)
 
template<>
std::string fixedSizeFormat (const std::string &s, const uint n)
 
template<typename T >
std::string vToString (const T &x)
 Convert something that can iterate into a string...
 
template<typename T >
std::string vectorAsStringWithCommas (const std::vector< T > &v)
 Convert a vector of type T to a string-list with commas.
 
template<>
std::string vectorAsStringWithCommas (const std::vector< std::string > &v)
 Specialization for strings that sanitizes the contained strings.
 
template<typename T >
std::vector< T > uniquifyVector (const std::vector< T > &list)
 Return a vector containing only the unique elements of the input vector.
 
base_generator_type & rng_singleton (void)
 Suite-wide random number generator singleton.
 
uint randomSeedRNG (void)
 Randomly seeds the RNG.
 
GCoord boxFromRemarks (const Remarks &)
 Extract the Alan-style box-size from a PDB Remarks block.
 
bool remarksHasBox (const Remarks &)
 Checks to see if a Remarks block has an Alan-style box size in it.
 
AtomicGroup loadStructureWithCoords (const std::string &model, const std::string &cooords)
 Loads a structure and optional coordinates.
 
AtomicGroup loadStructureWithCoords (const std::string &model_name, const std::string &type, const std::string &coord_name)
 
std::vector< uint > assignTrajectoryFrames (const pTraj &traj, const std::string &frame_index_spec, uint skip=0, uint stride=1)
 Builds a list of trajectory indices (frame_index_spec supercedes skip)
 

Variables

const std::string version_string = "LOOS 4.1.0 BUILT ON 20240212 16:31:52"
 
const uint kilobytes = 1024
 
const uint megabytes = 1024 * kilobytes
 
const uint gigabytes = 1024 * megabytes
 

Detailed Description

Namespace for most things not already encapsulated within a class.

Function Documentation

◆ averageStructure() [1/3]

AtomicGroup loos::averageStructure ( const AtomicGroup & ,
const std::vector< XForm > & ,
pTraj & traj )

Compute the average structure using all frames in a trajectory.

This version only reads a frame at a time from the trajectory. The trajectory iterator will be left pointing to the end of the trajectory.

◆ averageStructure() [2/3]

AtomicGroup loos::averageStructure ( const AtomicGroup & ,
const std::vector< XForm > & ,
pTraj & traj,
const std::vector< uint > & indices )

Compute the average structure from a trajectory reading only certain frames.

Note that the trajectory is NOT stored in memory. Frames will be read as needed. The trajectory "iterator" will be left pointing to the next frame after the last frame index passed.

◆ averageStructure() [3/3]

AtomicGroup loos::averageStructure ( const std::vector< AtomicGroup > & ensemble,
const std::vector< XForm > & xforms )

Compute the average structure of a set of AtomicGroup objects.

Takes into consideration the passed set of transforms...

◆ boxFromRemarks()

GCoord loos::boxFromRemarks ( const Remarks & )

Extract the Alan-style box-size from a PDB Remarks block.

Returns a GCoord(99999.99, 99999.99, 99999.99) if there is no box info found in the remarks block.

◆ createSystem()

AtomicGroup loos::createSystem ( const std::string & filename)

Factory function for reading in structure files.

This function will try to determine the filetype for a structure file by examining the suffix of the file. It will return an AtomicGroup copy of the input structure.

Amber files are handled differently... If the filename ends with ".prmtop", then LOOS will check for a corresponding ".inpcrd" file. If one is found, then the coords will be loaded into the group. Otherwise, the prmtop will be loaded without coords and returned.

◆ createTrajectory()

pTraj loos::createTrajectory ( const std::string & ,
const AtomicGroup &  )

Factory function for reading in a trajectory file.

This function will try to determine the filetype for a trajectory by examining the suffix of the file. It will return a boost shared pointer to a new Trajectory object.

It is very important to understand that the object returned by this function must behave polymorphically. That's why it is wrapped in a boost shared pointer to the base class. Do not try to deference it and assign it to a Trajectory object...

◆ deduceAtomicNumberFromMass()

unsigned int loos::deduceAtomicNumberFromMass ( const double mass,
const double tolerance = 0.1 )

Deduce an atomic number from the mass.

Only the first 96 elements are included in LOOS. An atomic number of 0 is returned if the mass is not found within LOOS' table.

◆ invocationHeader()

std::string loos::invocationHeader ( int ,
char * [] )

Create an invocation header.

This is a string that can be embedded in output that records the invoking user, command-line, and a timestamp.

◆ iterativeAlignment()

boost::tuple< std::vector< XForm >, greal, int > loos::iterativeAlignment ( const AtomicGroup & model,
pTraj & traj,
const std::vector< uint > & frame_indices,
greal threshold = 1e-6,
int maxiter = 1000 )

Compute an iterative superposition by reading in frames from the Trajectory.

The iterativeAlignment() functions that take a trajectory as an argument do NOT cache frames of the trajectory internally. This means that the trajectory will be read as many times as is necessary for the alignment to converge. In practice, the OS-specific caching will likely result in decent performance. If speed is essential, then consider using the iterativeAlignment() version that takes a std::vector<AtomicGroup>& as argument instead.

◆ operator*()

template<class T >
Coord< T > loos::operator* ( const Matrix44< T > & M,
const Coord< T > & v )

Matrix-vector multiply This has to be a friend outside the class for GCC to be happy...

◆ operator<<() [1/2]

std::ostream & loos::operator<< ( std::ostream & os,
const Gromacs & g )

Output the group as a GRO...

Output as a GRO.

◆ operator<<() [2/2]

std::ostream & loos::operator<< ( std::ostream & os,
const PDB & p )

Output the group as a PDB...

Output as a PDB.

There are some formatting changes that occur when the group has a large number of atoms or resids. The most significant is when you have 100,000 or more, in which case you lose the altloc and chainid fields on output. However, the output PDB will load into pymol...

◆ parseIndexRange()

std::vector< uint > loos::parseIndexRange ( const std::string & input,
const uint maxsize )

Use boost::Spirit to parse a new-style range list.

Similar to parseRangeList(), but uses boost::spirit can requires a maximum index size. This allows you to specify ranges without knowing the endpoint, e.g. 10:2:, means start at 10, skipping 2, until the value of maxsize. This is useful for trajectories, otherwise you would have to know how big they are before using a range to specify a skip and stride.

◆ parseRange()

template<typename T >
std::vector< T > loos::parseRange ( const std::string & text,
const T endpoint = 0 )

Parse an Octave/Matlab-style range.

The Octave format is one of the following:

  • value
  • start:stop
  • start:step:stop

The range is inclusive for both ends. Internally, a vector of T's is created. There is no checking to make sure that the vector doesn't completely fill up memory.

As with matlab/octave, to count down, the step must be negative and start > stop.

◆ parseRangeList()

template<typename T >
std::vector< T > loos::parseRangeList ( const std::string & text,
const T endpoint = 0 )

Parses a comma-separated list of Octave-style ranges.

This function breaks apart a string at the commas and passes each substring to parseRange. The union of all of the vectors returned by parseRange is then sorted in ascending order and duplicate values are removed. This vector is then returned to the caller.

◆ randomSeedRNG()

uint loos::randomSeedRNG ( void )

Randomly seeds the RNG.

Currently uses time(3) to seed the RNG obtained from the singleton... Returns the seed used.

◆ readTable()

template<typename T >
std::vector< std::vector< T > > loos::readTable ( LineReader & reader)

Read in a table of items using a LineReader object.

This is distinct from the Matrix class reader which requires a specific format. Here, comments may be freely interspersed with the data and rows may contain different numbers of columns

◆ readVector()

template<typename T >
std::vector< T > loos::readVector ( LineReader & reader)

Read a list of items using a LineReader object.

The readVector<T>() family of functions allows you to read in a list of items (numbers) from a stream or a file. When used with a LineReader object, you have control over how blank lines and comments are handled. When used with either an istream or a string, the default behavior is to skip blank lines and comments will begin with the '#' character and are stripped.

◆ rng_singleton()

base_generator_type & loos::rng_singleton ( void )

Suite-wide random number generator singleton.

 LOOS makes no assumptions about how the random number generator
 gets seeded.  It is up to the tool-writer to seed it with a known
 value,
Simple matrix template class using policy classes to determine behavior.
Definition MatrixImpl.hpp:148
base_generator_type & rng_singleton(void)
Suite-wide random number generator singleton.
Definition utils_random.cpp:27

or call randomSeedRNG() to randomly seed the random number generator...

◆ selectAtoms()

AtomicGroup loos::selectAtoms ( const AtomicGroup & source,
const std::string selection )

Applies a string-based selection to an atomic group...

This routine parses the passed string, turning it into a selector and applies it to source. If there is an exception in the parsing, this is repackaged into a more sensible error message (including the string that generated the error). No other exceptions are caught.

We're also assuming that you're always wanting to select some atoms, so lack of selection constitutes an error and an exception is thrown. Note that in both the case of a parse error and null-selection, a runtime_error exception is thrown so the catcher cannot disambiguate between the two.

◆ sortedIndex()

template<typename T , class SortPredicate >
std::vector< uint > loos::sortedIndex ( const T & A)

Sort a container using the given sort policy, returning the indices that permutes the container into the sorted order. The container to be sorted must support T::size() and T::operator[]. What is returned is a vector of unsigned ints that represent the index into the container when it is sorted.

◆ svd()

boost::tuple< RealMatrix, RealMatrix, RealMatrix > loos::svd ( std::vector< AtomicGroup > & ensemble,
const bool align = true )

Compute the SVD of an ensemble with optional alignment (note RSVs returned are transposed)

Returns the U, S, and V' of the SVD of the passed ensemble. If align is true, then the ensemble is iteratively aligned prior to computing the SVD.

◆ swab()

template<typename T >
T loos::swab ( const T & datum)

Returns a byte-swapped copy of an arbitrary type.

Only valid for simple types (i.e. int, float, double)

◆ writeAsciiMatrix() [1/3]

template<class T , class P , template< typename > class S>
void loos::writeAsciiMatrix ( const std::string & fname,
const Math::Matrix< T, P, S > & M,
const std::string & meta,
const bool trans = false )

Write an entire matrix to a file.

Principal front-end to writing out ASCII matrices. Given a filename, matrix, metadata, and optional flag for whether or not the matrix is transposed interally, will write the entire matrix using the default C++ type-to-string conversion.

◆ writeAsciiMatrix() [2/3]

template<class T , class P , template< typename > class S>
std::ostream & loos::writeAsciiMatrix ( std::ostream & os,
const Math::Matrix< T, P, S > & M,
const std::string & meta,
const bool trans = false )

Write an entire matrix to a stream.

 Primary front-end to writing an entire ASCII matrix to a stream.
 Will use the default C++ type-to-string conversion for the matrix
 data type (e.g. float).
 Example:
writeAsciiMatrix(cout, M, "Eigenvector matrix");
std::ostream & writeAsciiMatrix(std::ostream &os, const Math::Matrix< T, P, S > &M, const std::string &meta, const Math::Range &start, const Math::Range &end, const bool trans=false, F fmt=F())
Write a submatrix to a stream.
Definition MatrixWrite.hpp:140

◆ writeAsciiMatrix() [3/3]

template<class T , class P , template< typename > class S, class F >
std::ostream & loos::writeAsciiMatrix ( std::ostream & os,
const Math::Matrix< T, P, S > & M,
const std::string & meta,
const Math::Range & start,
const Math::Range & end,
const bool trans = false,
F fmt = F() )

Write a submatrix to a stream.

This family of functions write a matrix in ASCII format suitable for loading into Octave/Matlab or gnuplot. The meta information is written as part of the comment at the start of the file. The Math::Range start and end are just pairs that give an (j,i) starting and ending point within the matrix to write. Note that these arguments are not always honored (such as with a triangular matrix). The trans flag causes the output matrix to be the transpose of the stored matrix. The fmt arg is a functor that is expected for format each element of the matrix as a string. You can use this to adjust the precision of the output or delimit it, etc. The default is to use whatever the default operator<<() would be for type T.