n2p2 - A neural network potential package
Loading...
Searching...
No Matches
InterfaceLammps.h
Go to the documentation of this file.
1// n2p2 - A neural network potential package
2// Copyright (C) 2018 Andreas Singraber (University of Vienna)
3//
4// This program is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published by
6// the Free Software Foundation, either version 3 of the License, or
7// (at your option) any later version.
8//
9// This program is distributed in the hope that it will be useful,
10// but WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12// GNU General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with this program. If not, see <https://www.gnu.org/licenses/>.
16
17#ifndef INTERFACELAMMPS_H
18#define INTERFACELAMMPS_H
19
20#include "Mode.h"
21#include "Structure.h"
22#include <map> // std::map
23#include <cstddef> // std::size_t
24#include <cstdint> // int64_t
25#include <string> // std::string
26#include <vector> // std::vector
27
28namespace nnp
29{
30
31class InterfaceLammps : public Mode
32{
33public:
35
54 void initialize(char const* const& directory,
55 char const* const& emap,
56 bool showew,
57 bool resetew,
58 int showewsum,
59 int maxew,
60 double cflength,
61 double cfenergy,
62 double lammpsCutoff,
63 int lammpsNtypes,
64 int myRank);
69 void setGlobalStructureStatus(bool const status);
77 void setLocalAtoms(int numAtomsLocal,
78 int const* const atomType);
83 void setLocalAtomPositions(double const* const* const atomPos);
88 void setLocalTags(int const* const atomTag);
93 void setLocalTags(int64_t const* const atomTag);
102 void setBoxVectors(double const* boxlo,
103 double const* boxhi,
104 double const xy,
105 double const xz,
106 double const yz);
111 void allocateNeighborlists(int const* const numneigh);
125 void addNeighbor(int i,
126 int j,
127 int64_t tag,
128 int type,
129 double dx,
130 double dy,
131 double dz,
132 double d2);
140 void process();
144 void processDevelop();
149 double getEnergy() const;
158 double getAtomicEnergy(int index) const;
164 void addElectrostaticEnergy(double energy);
169 void getForces(double* const* const& atomF) const;
175 void getForcesDevelop(double* const* const& atomF) const;
180 void getForcesChi(double const* const& lambda,
181 double* const* const& atomF) const;
187 void getCharges(double* const& atomQ) const;
192 bool isInitialized() const;
197 double getMaxCutoffRadius() const;
202 double getEwaldPrec() const;
214 long getEWBufferSize() const;
220 void fillEWBuffer(char* const& buf, int bs) const;
226 void extractEWBuffer(char const* const& buf, int bs);
235 void addCharge(int index, double Q);
243 void getQEqParams(double* const& atomChi,
244 double* const& atomJ,
245 double* const& sigmaSqrtPi,
246 double *const *const& gammaSqrt2,
247 double& qRef) const;
253 void getdEdQ(double* const& dEtotdQ) const;
258 void getScreeningInfo(double* const& rScreen) const;
266 void getdChidxyz(int tag,
267 double* const& dChidx,
268 double* const& dChidy,
269 double* const& dChidz) const;
272 void setElecDone();
278 void writeToFile(std::string const fileName,
279 bool const append);
285 void add3DVecToArray(double *const & arr, Vec3D const& v) const;
286
287protected:
295 bool showew;
301 int maxew;
303 double cflength;
305 double cfenergy;
307 std::string emap;
309 std::vector<size_t> indexMap;
311 std::map<int, bool> ignoreType;
313 std::map<int, std::size_t> mapTypeToElement;
315 std::map<std::size_t, int> mapElementToType;
320};
321
323// Inlined function definitions //
325
327{
328 return initialized;
329}
330
331}
332
333#endif
Structure structure
Structure containing local atoms.
long getEWBufferSize() const
Calculate buffer size for extrapolation warning communication.
bool initialized
Initialization state.
void process()
Calculate symmetry functions, atomic neural networks and sum of local energy contributions.
void extractEWBuffer(char const *const &buf, int bs)
Extract given buffer to symmetry function statistics class.
void addCharge(int index, double Q)
Read atomic charges from LAMMPS into n2p2.
std::map< int, bool > ignoreType
True if atoms of this LAMMPS type will be ignored.
void setLocalTags(int const *const atomTag)
Set atom tags (int version, -DLAMMPS_SMALLBIG).
void allocateNeighborlists(int const *const numneigh)
Allocate neighbor lists.
void setElecDone()
Set isElecDone true after running the first NN in 4G-HDNNPs.
std::map< int, std::size_t > mapTypeToElement
Map from LAMMPS type to n2p2 element index.
void addElectrostaticEnergy(double energy)
Adds electrostatic energy contribution to the total structure energy.
void add3DVecToArray(double *const &arr, Vec3D const &v) const
Add a Vec3D vector to a 3D array in place.
void processDevelop()
Calculate symmetry functions, atomic neural networks and sum of local energy contributions (developme...
void getQEqParams(double *const &atomChi, double *const &atomJ, double *const &sigmaSqrtPi, double *const *const &gammaSqrt2, double &qRef) const
Write QEq arrays from n2p2 to LAMMPS.
double getMaxCutoffRadiusOverall()
Get largest cutoff including structure specific cutoff and screening cutoff.
bool hasGlobalStructure
Whether n2p2 knows about the global structure or only a local part.
bool isElecDone
True if first NN is calculated.
void setLocalAtomPositions(double const *const *const atomPos)
Set absolute atom positions from LAMMPS (nnp/develop only).
double getEnergy() const
Return sum of local energy contributions.
void fillEWBuffer(char *const &buf, int bs) const
Fill provided buffer with extrapolation warning entries.
double getAtomicEnergy(int index) const
Return energy contribution of one atom.
void getForces(double *const *const &atomF) const
Calculate forces and add to LAMMPS atomic force arrays.
bool resetew
Corresponds to LAMMPS resetew keyword.
bool isInitialized() const
Check if this interface is correctly initialized.
void writeExtrapolationWarnings()
Write extrapolation warnings to log.
void setBoxVectors(double const *boxlo, double const *boxhi, double const xy, double const xz, double const yz)
Set box vectors of structure stored in LAMMPS (nnp/develop only).
std::vector< size_t > indexMap
Map from LAMMPS index to n2p2 atom index.
bool getGlobalStructureStatus()
Check if n2p2 knows about global structure.
double cfenergy
Corresponds to LAMMPS cfenergy keyword.
void getScreeningInfo(double *const &rScreen) const
Read screening function information from n2p2 into LAMMPS.
void setLocalAtoms(int numAtomsLocal, int const *const atomType)
(Re)set structure to contain only local LAMMPS atoms.
bool showew
Corresponds to LAMMPS showew keyword.
void writeToFile(std::string const fileName, bool const append)
Write current structure to file in units used in training data.
void getCharges(double *const &atomQ) const
Transfer charges (in units of e) to LAMMPS atomic charge vector.
int showewsum
Corresponds to LAMMPS showewsum keyword.
void setGlobalStructureStatus(bool const status)
Specify whether n2p2 knows about global structure or only local structure.
void finalizeNeighborList()
Sorts neighbor list and creates cutoff map if necessary.
double getEwaldPrec() const
Get Ewald precision parameter.
void getForcesDevelop(double *const *const &atomF) const
Calculate forces and add to LAMMPS atomic force arrays (development version for "hdnnp/develop" pair ...
void clearExtrapolationWarnings()
Clear extrapolation warnings storage.
void getForcesChi(double const *const &lambda, double *const *const &atomF) const
Calculate chi-term for forces and add to LAMMPS atomic force arrays.
std::map< std::size_t, int > mapElementToType
Map from n2p2 element index to LAMMPS type.
void getdEdQ(double *const &dEtotdQ) const
Write the derivative of total energy with respect to atomic charges from n2p2 into LAMMPS.
int myRank
Process rank.
double getMaxCutoffRadius() const
Get largest cutoff of symmetry functions.
void getdChidxyz(int tag, double *const &dChidx, double *const &dChidy, double *const &dChidz) const
Transfer spatial derivatives of atomic electronegativities.
int maxew
Corresponds to LAMMPS maxew keyword.
void addNeighbor(int i, int j, int64_t tag, int type, double dx, double dy, double dz, double d2)
Add one neighbor to atom (int64_t version, -DLAMMPS_BIGBIG).
double cflength
Corresponds to LAMMPS cflength keyword.
std::string emap
Corresponds to LAMMPS map keyword.
void initialize()
Write welcome message with version information.
Definition Mode.cpp:56
Definition Atom.h:29
Storage for one atomic configuration.
Definition Structure.h:39
Vector in 3 dimensional real space.
Definition Vec3D.h:30