atomicrex  0.1
An advanced atomistic model building tool
Public Types | Public Member Functions | Protected Attributes | List of all members
atomicrex::AtomicStructure Class Reference

Base class for maintaining structures. More...

#include <AtomicStructure.h>

Inheritance diagram for atomicrex::AtomicStructure:
atomicrex::FitObject atomicrex::BinaryCubicLatticeStructure atomicrex::BinaryHexaRhomboTetraLatticeStructure atomicrex::DimerStructure atomicrex::SuperCellStructure atomicrex::UnaryCubicLatticeStructure atomicrex::UnaryHexaRhomboTetraLatticeStructure atomicrex::UserStructure atomicrex::B1LatticeStructure atomicrex::B2LatticeStructure atomicrex::B3LatticeStructure atomicrex::C15LatticeStructure atomicrex::C1LatticeStructure atomicrex::D8aLatticeStructure atomicrex::L12LatticeStructure atomicrex::B4LatticeStructure atomicrex::BhLatticeStructure atomicrex::D2dLatticeStructure atomicrex::L10LatticeStructure atomicrex::Ni17Th2LatticeStructure atomicrex::Th2Zn17LatticeStructure atomicrex::PointDefectStructure atomicrex::BCCLatticeStructure atomicrex::DIALatticeStructure atomicrex::FCCLatticeStructure atomicrex::SCLatticeStructure atomicrex::betaSnLatticeStructure atomicrex::DHCPLatticeStructure atomicrex::HCPLatticeStructure atomicrex::OMGLatticeStructure

Public Types

enum  DirtyFlags { NEIGHBOR_LISTS = (1<<0), ATOM_POSITIONS = (1<<1), STRUCTURE_DOF = (1<<4) }
 Flags that indicate which parts of the structure need to be updated prior to the next energy calculation. More...
 
enum  OutputFormat { POSCAR, LAMMPS }
 Possible formats for writing the structure to file.
 

Public Member Functions

 AtomicStructure (const FPString &id, FitJob *job)
 
void setupSimulationCell (const Matrix3 &cellVectors, const Point3 &cellOrigin=Point3::Origin(), const std::array< bool, 3 > &pbc=std::array< bool, 3 >{{true, true, true}})
 Set up the simulation cell. More...
 
void deformSimulationCell (const Matrix3 &deformation)
 Applies an affine transformation to the cell and the atoms. More...
 
const Matrix3 & simulationCell () const
 Returns global simulation cell matrix. More...
 
const Point3 & simulationCellOrigin () const
 The origin point of the global simulation cell in world coordinates. More...
 
const Matrix3 & reciprocalSimulationCell () const
 The inverse of the simulation cell matrix used to transform absolute coordinates to reduced coordinates. More...
 
bool hasPBC (size_t dimension) const
 Returns whether periodic boundary conditions are enabled in the given spatial dimension. More...
 
const std::array< bool, 3 > & pbc () const
 Returns the periodic boundary condition flags. More...
 
Vector3I periodicImage (const Point3 &p) const
 Determines the periodic image of the cell the given point is located in.
 
Point3 wrapPoint (Point3 p) const
 Wraps a point to be inside the simulation cell if periodic boundary conditions are enabled.
 
Point3 wrapReducedPoint (Point3 p) const
 Wraps a point given in reduced coordinated to be inside the simulation cell.
 
Point3 reducedToAbsolute (const Point3 &reducedPoint) const
 Converts a point given in reduced cell coordinates to a point in absolute coordinates.
 
Point3 absoluteToReduced (const Point3 &worldPoint) const
 Converts a point given in absolute coordinates to a point in reduced cell coordinates.
 
Vector3 reducedToAbsolute (const Vector3 &reducedVec) const
 Converts a vector given in reduced cell coordinates to a vector in absolute coordinates.
 
Vector3 absoluteToReduced (const Vector3 &worldVec) const
 Converts a vector given in absolute coordinates to a point in vector cell coordinates.
 
void setAtomCount (int numLocalAtoms)
 Resizes the atoms array. More...
 
int numLocalAtoms () const
 Returns the number of (real) atoms in the structure cell.
 
int numGhostAtoms () const
 Returns the number of ghost atoms in the structure cell.
 
int numAtoms () const
 Returns the total number of atoms in the structure cell including real and ghost atoms.
 
const std::vector< Point3 > & atomPositions () const
 Returns the array of positions of all atoms (including ghosts) in the structure (const version).
 
std::vector< Point3 > & atomPositions ()
 Returns the array of positions of all atoms (including ghosts) in the structure cell (non-const version).
 
void setAtomPositions (const std::vector< Point3 > &newPositions)
 Sets the array of positions of all atoms in the structure cell. More...
 
const std::vector< Point3 > & atomDisplacements () const
 Returns the array of displacements of all atoms in the structure cell wrt a reference structure (const version).
 
std::vector< Point3 > & atomDisplacements ()
 Returns the array of displacements of all atoms in the structure cell wrt a reference structure (non-const version).
 
const std::vector< Point3 > & atomInitialPositions () const
 Returns the initial positions of the all real atoms they had at the beginning of the job.
 
int atomType (int i) const
 Returns the type of a i-th atom.
 
const std::vector< int > & atomTypes () const
 Returns the array of types of all atoms in the cell (const version).
 
std::vector< int > & atomTypes ()
 Returns the array of types of all atoms in the cell (non-const version).
 
const std::vector< int > & atomTags () const
 Returns the array of IDs of the atoms in the cell (const version).
 
std::vector< int > & atomTags ()
 Returns the array of IDs of the atoms in the cell (non-const version).
 
const std::vector< Vector3 > & atomForces () const
 Returns the array of force vectors (const version).
 
std::vector< Vector3 > & atomForces ()
 Returns the array of force vectors (non-const version).
 
AtomVectorPropertyatomForcesProperty ()
 Returns the property for the atomic force vectors.
 
const std::vector< std::pair< int, int > > & forwardMapping () const
 Returns an array of integer pairs that map the real atoms to the ghost atoms.
 
const std::vector< int > & reverseMapping () const
 Returns an array of integers that maps the ghost atoms to the real atoms.
 
template<typename T >
T * perAtomData ()
 Returns a pointer to the memory buffer that stores temporary per-atom data used by the potential routines.
 
NeighborListneighborList (const Potential *potential)
 
void setDirty (DirtyFlags flags)
 Marks parts of the structure that must be updated before the next energy calculation.
 
bool isDirty (DirtyFlags parts) const
 Tests whether certain parts of the structure must be updated before the next energy calculation.
 
void clearDirty (DirtyFlags parts)
 Resets the dirty flags.
 
virtual void dofValueChanged (DegreeOfFreedom &dof) override
 This callback function is called by the DOFs of the structure each time when their values changes. More...
 
virtual void updateStructure ()
 Updates the structure, creates ghost atoms and builds neighbor lists.
 
virtual bool computeProperties (bool isFitting) override
 Computes all enabled properties of the structures. More...
 
double computeEnergy (bool computeForces, bool isFitting, bool suppressRelaxation=false)
 Computes the total energy and optionally the forces for this structure. More...
 
bool relax (bool isFitting)
 Relaxes the structural degrees of freedom such that the total energy is minimized. More...
 
void writeToDumpFile (const FPString &filename, bool includeGhostAtoms=false) const
 Exports the structure to a LAMMPS dump file. More...
 
void writeToPoscarFile (const FPString &filename, bool includeGhostAtoms=false) const
 Exports the structure to a POSCAR file. More...
 
void writeToFile ()
 Exports the structure to file if a filename has been provided.
 
double totalEnergy () const
 Returns the total potential energy of this structure after computeEnergy() has been called.
 
const std::array< double, 6 > & virial () const
 Returns the virial tensor computed by the force routine.
 
std::array< double, 6 > & virial ()
 Returns a reference to the virial tensor. More...
 
double pressureTensor (int voigtIndex) const
 Returns a component of the pressure tensor. More...
 
double pressure () const
 Returns the hydrostatic part of the pressure tensor. More...
 
virtual void parse (XML::Element structureElement) override
 Parses the structure-specific parameters in the XML element in the job file.
 
- Public Member Functions inherited from atomicrex::FitObject
virtual ~FitObject ()=default
 Virtual destructor.
 
double relativeWeight () const
 Returns the relative fit weight assigned to this object.
 
void setRelativeWeight (double weight)
 Assigns a relative fit weight to this object.
 
virtual void assignAbsoluteWeights (double absoluteWeight)
 Recursively assigns absolute weights to the properties of this object and its sub-objects.
 
const std::vector< FitProperty * > & properties () const
 Returns a list of fitting properties of this object.
 
void listAllProperties (std::vector< FitProperty *> &list) const
 Builds a list of properties of this object and all its sub-objects.
 
FitPropertypropertyById (const FPString &id) const
 Returns the property with the given ID.
 
const std::vector< DegreeOfFreedom * > & DOF () const
 Returns a list of degrees of freedom of this object.
 
void listAllDOF (std::vector< DegreeOfFreedom *> &list) const
 Builds a list of degrees of freedom of this object and all its sub-objects.
 
DegreeOfFreedomDOFById (const FPString &id, const FPString &tag=FPString()) const
 Returns the degree of freedom with the given ID (and tag).
 
const std::vector< FitObject * > & fitObjects () const
 Returns the list of FitObjects which are part of this group.
 
void registerSubObject (FitObject *subobject, bool deleteOnShutdown=false)
 Registers a sub-object.
 
virtual void print (MsgLogger &stream)
 Outputs the name of the object.
 
bool fitEnabled () const
 Returns whether this object and it's children are included in the fit.
 
void setFitEnabled (bool enable)
 Sets whether this object and it's children are included in the fit.
 
bool outputEnabled () const
 
void setOutputEnabled (bool enable)
 
const FPStringid () const
 Returns the identifier of this object instance.
 
void setId (const FPString &id)
 Sets the main identification tag.
 
const FPStringtag () const
 Returns the assigned tag string.
 
void setTag (const FPString &tag)
 Sets the complementary identification tag.
 
FitJobjob () const
 Returns a pointer to the job to which this object belongs.
 
FitObjectparent () const
 Returns the parent of this object in the hierarchy.
 

Protected Attributes

ScalarFitProperty _totalEnergyProperty
 Output of the energy calculation that can be used for fitting.
 
ScalarFitProperty _atomicEnergyProperty
 Output of the per atom energy calculation that can be used for fitting.
 
ScalarFitProperty _totalVolumeProperty
 The total volume of the structure cell.
 
ScalarFitProperty _atomicVolumeProperty
 The average volume per atom in the structure cell.
 
ScalarFitProperty _pressureTensorProperty [6]
 The components of the cell's pressure tensor.
 
ScalarFitProperty _pressureProperty
 One third of the trace of the cell's pressure tensor.
 
ScalarFitProperty _bulkModulusProperty
 The bulk modulus of the lattice.
 
ScalarFitProperty _elasticConstantProperties [21]
 The elastic constants.
 
AtomVectorProperty _atomForcesProperty
 The force vector.
 
AtomCoordinatesDOF _atomCoordinatesDOF
 This DOF encapsulates the atomic degrees of freedom.
 
- Protected Attributes inherited from atomicrex::FitObject
bool _fitEnabled = true
 Controls whether this object and it's children are included in the fit.
 
bool _outputEnabled = true
 
FitJob_job = nullptr
 Pointer to the job this object belongs to.
 
FPString _id
 The identifier string of this object instance.
 
double _relativeWeight = 1.0
 The relative fit weight assigned to this object.
 

Additional Inherited Members

- Protected Member Functions inherited from atomicrex::FitObject
 FitObject ()
 Default Constructor.
 
 FitObject (const FPString &id, FitJob *job, const FPString &tag=FPString())
 Constructor.
 
void registerProperty (FitProperty *prop, bool deleteOnShutdown=false)
 Registers a property of this object.
 
void registerDOF (DegreeOfFreedom *dof)
 Registers a DOF of this object.
 

Detailed Description

Base class for maintaining structures.

This class is used to store structure specific information including atomic coordinates, neighbor lists, total energy, volume, pressure etc.

Member Enumeration Documentation

◆ DirtyFlags

Flags that indicate which parts of the structure need to be updated prior to the next energy calculation.

Enumerator
NEIGHBOR_LISTS 

neighbor lists require updating

ATOM_POSITIONS 

ghost atom positions require updating

STRUCTURE_DOF 

degree(s) of freedom requires updating

Member Function Documentation

◆ computeEnergy()

double atomicrex::AtomicStructure::computeEnergy ( bool  computeForces,
bool  isFitting,
bool  suppressRelaxation = false 
)

Computes the total energy and optionally the forces for this structure.

Computes the total energy and optionally the forces of this structure.

Parameters
computeForcesTurns on the computation of forces.
isFittingIndicates that the structure is included in the fitting data set.
suppressRelaxationIf True structure relaxation will be suppressed.
Returns
A scalar representing the energy of the structure. If the forces are computed they are stored in the structure object.
See also
atomForcesProperty() for retrieving the atomic forces

By default, performs relaxation of atomic positions if enabled by the user.

Parameters
computeForcesif True the forces will be evaluated as well
isFittingdetermines whether this DOF should be relaxed (usually the desired behavior depends on the program phase)
suppressRelaxationif True structural relaxation will be suppressed

◆ computeProperties()

bool atomicrex::AtomicStructure::computeProperties ( bool  isFitting)
overridevirtual

Computes all enabled properties of the structures.

Parameters
isFittingIndicates that the structure is included in the fitting data set.

Reimplemented from atomicrex::FitObject.

◆ deformSimulationCell()

void atomicrex::AtomicStructure::deformSimulationCell ( const Matrix3 &  deformation)

Applies an affine transformation to the cell and the atoms.

Parameters
deformationdefined 3x3 deformation matrix that is applied to the simulation cell

◆ dofValueChanged()

virtual void atomicrex::AtomicStructure::dofValueChanged ( DegreeOfFreedom dof)
inlineoverridevirtual

This callback function is called by the DOFs of the structure each time when their values changes.

Warning
This method should not be called by the user.

Reimplemented from atomicrex::FitObject.

◆ hasPBC()

bool atomicrex::AtomicStructure::hasPBC ( size_t  dimension) const
inline

Returns whether periodic boundary conditions are enabled in the given spatial dimension.

Parameters
dimensionindex to spatial dimension (0=x, 1=y, z=3)
Returns
True/False if boundary condition is enabled/disabled.

◆ pbc()

const std::array<bool,3>& atomicrex::AtomicStructure::pbc ( ) const
inline

Returns the periodic boundary condition flags.

Returns
3-dimensional vector of True/False values.

◆ pressure()

double atomicrex::AtomicStructure::pressure ( ) const
inline

Returns the hydrostatic part of the pressure tensor.

The function thus returns one third of the trace of the pressure tensor.

◆ pressureTensor()

double atomicrex::AtomicStructure::pressureTensor ( int  voigtIndex) const
inline

Returns a component of the pressure tensor.

Parameters
voigtIndexIndex of the pressure tensor component in Voigt notation

◆ reciprocalSimulationCell()

const Matrix3& atomicrex::AtomicStructure::reciprocalSimulationCell ( ) const
inline

The inverse of the simulation cell matrix used to transform absolute coordinates to reduced coordinates.

Returns
3x3 matrix representing the reciprocal cell metric

◆ relax()

bool atomicrex::AtomicStructure::relax ( bool  isFitting)

Relaxes the structural degrees of freedom such that the total energy is minimized.

Parameters
isFittingIndicates that the structure is included in the fitting data set.

◆ setAtomCount()

void atomicrex::AtomicStructure::setAtomCount ( int  numLocalAtoms)

Resizes the atoms array.

Parameters
numLocalAtomsnew number of (local) atoms.

◆ setAtomPositions()

void atomicrex::AtomicStructure::setAtomPositions ( const std::vector< Point3 > &  newPositions)
inline

Sets the array of positions of all atoms in the structure cell.

Parameters
newPositionsvector of 3-dimensional vectors ("points") representing the atomic coordinates. Size must match the number of local atoms.

◆ setupSimulationCell()

void atomicrex::AtomicStructure::setupSimulationCell ( const Matrix3 &  cellVectors,
const Point3 &  cellOrigin = Point3::Origin(),
const std::array< bool, 3 > &  pbc = std::array<bool,3>{{true,true,true}} 
)

Set up the simulation cell.

Parameters
cellVectorscell metric (3x3 matrix)
cellOriginorigin of cell (3-dim vector)
pbcperiodic boundary conditions in each direction

◆ simulationCell()

const Matrix3& atomicrex::AtomicStructure::simulationCell ( ) const
inline

Returns global simulation cell matrix.

Returns
3x3 matrix representing the simulation cell metric

◆ simulationCellOrigin()

const Point3& atomicrex::AtomicStructure::simulationCellOrigin ( ) const
inline

The origin point of the global simulation cell in world coordinates.

Returns
3-dimensional vector representing the origin of the simulation cell.

◆ virial()

std::array<double,6>& atomicrex::AtomicStructure::virial ( )
inline

Returns a reference to the virial tensor.

The virial tensor is returned as a six-dimensional vector. The indices of the vector follow the Voigt notation.

◆ writeToDumpFile()

void atomicrex::AtomicStructure::writeToDumpFile ( const FPString filename,
bool  includeGhostAtoms = false 
) const

Exports the structure to a LAMMPS dump file.

Parameters
filenameName of the output file
includeGhostAtomsIf True ghost atoms used during computation of the energy will be included in the output

◆ writeToPoscarFile()

void atomicrex::AtomicStructure::writeToPoscarFile ( const FPString filename,
bool  includeGhostAtoms = false 
) const

Exports the structure to a POSCAR file.

Parameters
filenameName of the output file
includeGhostAtomsIf True ghost atoms used during computation of the energy will be included in the output

The documentation for this class was generated from the following files: