n2p2 - A neural network potential package
Loading...
Searching...
No Matches
nnp::InterfaceLammps Class Reference

#include <InterfaceLammps.h>

Inheritance diagram for nnp::InterfaceLammps:
Collaboration diagram for nnp::InterfaceLammps:

Public Member Functions

 InterfaceLammps ()
 
void initialize (char const *const &directory, char const *const &emap, bool showew, bool resetew, int showewsum, int maxew, double cflength, double cfenergy, double lammpsCutoff, int lammpsNtypes, int myRank)
 Initialize the LAMMPS interface.
 
void setGlobalStructureStatus (bool const status)
 Specify whether n2p2 knows about global structure or only local structure.
 
bool getGlobalStructureStatus ()
 Check if n2p2 knows about global structure.
 
void setLocalAtoms (int numAtomsLocal, int const *const atomType)
 (Re)set structure to contain only local LAMMPS atoms.
 
void setLocalAtomPositions (double const *const *const atomPos)
 Set absolute atom positions from LAMMPS (nnp/develop only).
 
void setLocalTags (int const *const atomTag)
 Set atom tags (int version, -DLAMMPS_SMALLBIG).
 
void setLocalTags (int64_t const *const atomTag)
 Set atom tags (int64_t version, -DLAMMPS_BIGBIG).
 
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).
 
void allocateNeighborlists (int const *const numneigh)
 Allocate neighbor lists.
 
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).
 
void finalizeNeighborList ()
 Sorts neighbor list and creates cutoff map if necessary.
 
void process ()
 Calculate symmetry functions, atomic neural networks and sum of local energy contributions.
 
void processDevelop ()
 Calculate symmetry functions, atomic neural networks and sum of local energy contributions (development version for "hdnnp/develop" pair style).
 
double getEnergy () const
 Return sum of local energy contributions.
 
double getAtomicEnergy (int index) const
 Return energy contribution of one atom.
 
void addElectrostaticEnergy (double energy)
 Adds electrostatic energy contribution to the total structure energy.
 
void getForces (double *const *const &atomF) const
 Calculate forces and add to LAMMPS atomic force arrays.
 
void getForcesDevelop (double *const *const &atomF) const
 Calculate forces and add to LAMMPS atomic force arrays (development version for "hdnnp/develop" pair style).
 
void getForcesChi (double const *const &lambda, double *const *const &atomF) const
 Calculate chi-term for forces and add to LAMMPS atomic force arrays.
 
void getCharges (double *const &atomQ) const
 Transfer charges (in units of e) to LAMMPS atomic charge vector.
 
bool isInitialized () const
 Check if this interface is correctly initialized.
 
double getMaxCutoffRadius () const
 Get largest cutoff of symmetry functions.
 
double getEwaldPrec () const
 Get Ewald precision parameter.
 
double getMaxCutoffRadiusOverall ()
 Get largest cutoff including structure specific cutoff and screening cutoff.
 
long getEWBufferSize () const
 Calculate buffer size for extrapolation warning communication.
 
void fillEWBuffer (char *const &buf, int bs) const
 Fill provided buffer with extrapolation warning entries.
 
void extractEWBuffer (char const *const &buf, int bs)
 Extract given buffer to symmetry function statistics class.
 
void writeExtrapolationWarnings ()
 Write extrapolation warnings to log.
 
void clearExtrapolationWarnings ()
 Clear extrapolation warnings storage.
 
void addCharge (int index, double Q)
 Read atomic charges from LAMMPS into n2p2.
 
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.
 
void getdEdQ (double *const &dEtotdQ) const
 Write the derivative of total energy with respect to atomic charges from n2p2 into LAMMPS.
 
void getScreeningInfo (double *const &rScreen) const
 Read screening function information from n2p2 into LAMMPS.
 
void getdChidxyz (int tag, double *const &dChidx, double *const &dChidy, double *const &dChidz) const
 Transfer spatial derivatives of atomic electronegativities.
 
void setElecDone ()
 Set isElecDone true after running the first NN in 4G-HDNNPs.
 
void writeToFile (std::string const fileName, bool const append)
 Write current structure to file in units used in training data.
 
void add3DVecToArray (double *const &arr, Vec3D const &v) const
 Add a Vec3D vector to a 3D array in place.
 
- Public Member Functions inherited from nnp::Mode
 Mode ()
 
void initialize ()
 Write welcome message with version information.
 
void loadSettingsFile (std::string const &fileName="input.nn")
 Open settings file and load all keywords into memory.
 
void setupGeneric (std::string const &nnpDir="", bool skipNormalize=false, bool initialHardness=false)
 Combine multiple setup routines and provide a basic NNP setup.
 
void setupNormalization (bool standalone=true)
 Set up normalization.
 
virtual void setupElementMap ()
 Set up the element map.
 
virtual void setupElements ()
 Set up all Element instances.
 
void setupCutoff ()
 Set up cutoff function for all symmetry functions.
 
virtual void setupSymmetryFunctions ()
 Set up all symmetry functions.
 
void setupSymmetryFunctionScalingNone ()
 Set up "empy" symmetry function scaling.
 
virtual void setupSymmetryFunctionScaling (std::string const &fileName="scaling.data")
 Set up symmetry function scaling from file.
 
virtual void setupSymmetryFunctionGroups ()
 Set up symmetry function groups.
 
virtual void setupSymmetryFunctionCache (bool verbose=false)
 Set up symmetry function cache.
 
void setupSymmetryFunctionMemory (bool verbose=false)
 Extract required memory dimensions for symmetry function derivatives.
 
void setupSymmetryFunctionStatistics (bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
 Set up symmetry function statistics collection.
 
void setupCutoffMatrix ()
 Setup matrix storing all symmetry function cut-offs for each element.
 
virtual void setupNeuralNetwork ()
 Set up neural networks for all elements.
 
virtual void setupNeuralNetworkWeights (std::map< std::string, std::string > fileNameFormats=std::map< std::string, std::string >())
 Set up neural network weights from files with given name format.
 
virtual void setupNeuralNetworkWeights (std::string directoryPrefix, std::map< std::string, std::string > fileNameFormats=std::map< std::string, std::string >())
 Set up neural network weights from files with given name format.
 
virtual void setupElectrostatics (bool initialHardness=false, std::string directoryPrefix="", std::string fileNameFormat="hardness.%03zu.data")
 Set up electrostatics related stuff (hardness, screening, ...).
 
void calculateSymmetryFunctions (Structure &structure, bool const derivatives)
 Calculate all symmetry functions for all atoms in given structure.
 
void calculateSymmetryFunctionGroups (Structure &structure, bool const derivatives)
 Calculate all symmetry function groups for all atoms in given structure.
 
void calculateAtomicNeuralNetworks (Structure &structure, bool const derivatives, std::string id="")
 Calculate a single atomic neural network for a given atom and nn type.
 
void chargeEquilibration (Structure &structure, bool const derivativesElec)
 Perform global charge equilibration method.
 
void calculateEnergy (Structure &structure) const
 Calculate potential energy for a given structure.
 
void calculateCharge (Structure &structure) const
 Calculate total charge for a given structure.
 
void calculateForces (Structure &structure) const
 Calculate forces for all atoms in given structure.
 
void evaluateNNP (Structure &structure, bool useForces=true, bool useDEdG=true)
 Evaluate neural network potential (includes total energy, optionally forces and in some cases charges.
 
void addEnergyOffset (Structure &structure, bool ref=true)
 Add atomic energy offsets to reference energy.
 
void removeEnergyOffset (Structure &structure, bool ref=true)
 Remove atomic energy offsets from reference energy.
 
double getEnergyOffset (Structure const &structure) const
 Get atomic energy offset for given structure.
 
double getEnergyWithOffset (Structure const &structure, bool ref=true) const
 Add atomic energy offsets and return energy.
 
double normalized (std::string const &property, double value) const
 Apply normalization to given property.
 
double normalizedEnergy (Structure const &structure, bool ref=true) const
 Apply normalization to given energy of structure.
 
double physical (std::string const &property, double value) const
 Undo normalization for a given property.
 
double physicalEnergy (Structure const &structure, bool ref=true) const
 Undo normalization for a given energy of structure.
 
void convertToNormalizedUnits (Structure &structure) const
 Convert one structure to normalized units.
 
void convertToPhysicalUnits (Structure &structure) const
 Convert one structure to physical units.
 
void logEwaldCutoffs ()
 Logs Ewald params whenever they change.
 
std::size_t getNumExtrapolationWarnings () const
 Count total number of extrapolation warnings encountered for all elements and symmetry functions.
 
void resetExtrapolationWarnings ()
 Erase all extrapolation warnings and reset counters.
 
NNPType getNnpType () const
 Getter for Mode::nnpType.
 
double getMeanEnergy () const
 Getter for Mode::meanEnergy.
 
double getConvEnergy () const
 Getter for Mode::convEnergy.
 
double getConvLength () const
 Getter for Mode::convLength.
 
double getConvCharge () const
 Getter for Mode::convCharge.
 
double getEwaldPrecision () const
 Getter for Mode::ewaldSetup.precision.
 
double getEwaldMaxCharge () const
 Getter for Mode::ewaldSetup.maxCharge.
 
double getEwaldMaxSigma () const
 Getter for Mode::ewaldSetup.maxQsigma.
 
EWALDTruncMethod getEwaldTruncationMethod () const
 Getter for Mode::ewaldSetup.truncMethod.
 
KSPACESolver kspaceSolver () const
 Getter for Mode::kspaceSolver.
 
double getMaxCutoffRadius () const
 Getter for Mode::maxCutoffRadius.
 
std::size_t getNumElements () const
 Getter for Mode::numElements.
 
ScreeningFunction getScreeningFunction () const
 Getter for Mode::screeningFunction.
 
std::vector< std::size_t > getNumSymmetryFunctions () const
 Get number of symmetry functions per element.
 
bool useNormalization () const
 Check if normalization is enabled.
 
bool settingsKeywordExists (std::string const &keyword) const
 Check if keyword was found in settings file.
 
std::string settingsGetValue (std::string const &keyword) const
 Get value for given keyword in Settings instance.
 
std::vector< std::size_t > pruneSymmetryFunctionsRange (double threshold)
 Prune symmetry functions according to their range and write settings file.
 
std::vector< std::size_t > pruneSymmetryFunctionsSensitivity (double threshold, std::vector< std::vector< double > > sensitivity)
 Prune symmetry functions with sensitivity analysis data.
 
void writePrunedSettingsFile (std::vector< std::size_t > prune, std::string fileName="output.nn") const
 Copy settings file but comment out lines provided.
 
void writeSettingsFile (std::ofstream *const &file) const
 Write complete settings file.
 

Protected Attributes

int myRank
 Process rank.
 
bool initialized
 Initialization state.
 
bool hasGlobalStructure
 Whether n2p2 knows about the global structure or only a local part.
 
bool showew
 Corresponds to LAMMPS showew keyword.
 
bool resetew
 Corresponds to LAMMPS resetew keyword.
 
int showewsum
 Corresponds to LAMMPS showewsum keyword.
 
int maxew
 Corresponds to LAMMPS maxew keyword.
 
double cflength
 Corresponds to LAMMPS cflength keyword.
 
double cfenergy
 Corresponds to LAMMPS cfenergy keyword.
 
std::string emap
 Corresponds to LAMMPS map keyword.
 
std::vector< size_t > indexMap
 Map from LAMMPS index to n2p2 atom index.
 
std::map< int, bool > ignoreType
 True if atoms of this LAMMPS type will be ignored.
 
std::map< int, std::size_t > mapTypeToElement
 Map from LAMMPS type to n2p2 element index.
 
std::map< std::size_t, int > mapElementToType
 Map from n2p2 element index to LAMMPS type.
 
Structure structure
 Structure containing local atoms.
 
bool isElecDone
 True if first NN is calculated.
 
- Protected Attributes inherited from nnp::Mode
NNPType nnpType
 
bool normalize
 
bool checkExtrapolationWarnings
 
std::size_t numElements
 
std::vector< std::size_t > minNeighbors
 
std::vector< double > minCutoffRadius
 
double maxCutoffRadius
 
double cutoffAlpha
 
double meanEnergy
 
double convEnergy
 
double convLength
 
double convCharge
 
double fourPiEps
 
EwaldSetup ewaldSetup
 
KspaceGrid kspaceGrid
 
settings::Settings settings
 
SymFnc::ScalingType scalingType
 
CutoffFunction::CutoffType cutoffType
 
ScreeningFunction screeningFunction
 
std::vector< Elementelements
 
std::vector< std::string > nnk
 
std::map< std::string, NNSetupnns
 
std::vector< std::vector< double > > cutoffs
 Matrix storing all symmetry function cut-offs for all elements.
 
ErfcBuf erfcBuf
 

Additional Inherited Members

- Public Types inherited from nnp::Mode
enum class  NNPType { HDNNP_2G = 2 , HDNNP_4G = 4 , HDNNP_Q = 10 }
 
- Public Attributes inherited from nnp::Mode
ElementMap elementMap
 Global element map, populated by setupElementMap().
 
Log log
 Global log file.
 
- Protected Member Functions inherited from nnp::Mode
void readNeuralNetworkWeights (std::string const &id, std::string const &fileName)
 Read in weights for a specific type of neural network.
 

Detailed Description

Definition at line 31 of file InterfaceLammps.h.

Constructor & Destructor Documentation

◆ InterfaceLammps()

InterfaceLammps::InterfaceLammps ( )

Definition at line 38 of file InterfaceLammps.cpp.

38 : myRank (0 ),
39 initialized (false),
40 hasGlobalStructure (false),
41 showew (true ),
42 resetew (false),
43 showewsum (0 ),
44 maxew (0 ),
45 cflength (1.0 ),
46 cfenergy (1.0 ),
47 isElecDone (false)
48{
49}
bool initialized
Initialization state.
bool hasGlobalStructure
Whether n2p2 knows about the global structure or only a local part.
bool isElecDone
True if first NN is calculated.
bool resetew
Corresponds to LAMMPS resetew keyword.
double cfenergy
Corresponds to LAMMPS cfenergy keyword.
bool showew
Corresponds to LAMMPS showew keyword.
int showewsum
Corresponds to LAMMPS showewsum keyword.
int myRank
Process rank.
int maxew
Corresponds to LAMMPS maxew keyword.
double cflength
Corresponds to LAMMPS cflength keyword.

References cfenergy, cflength, hasGlobalStructure, initialized, isElecDone, maxew, myRank, resetew, showew, and showewsum.

Member Function Documentation

◆ initialize()

void InterfaceLammps::initialize ( char const *const & directory,
char const *const & emap,
bool showew,
bool resetew,
int showewsum,
int maxew,
double cflength,
double cfenergy,
double lammpsCutoff,
int lammpsNtypes,
int myRank )

Initialize the LAMMPS interface.

Parameters
[in]directoryDirectory containing NNP data files (weights, scaling, settings).
[in]emapElement mapping from LAMMPS to n2p2.
[in]showewIf detailed extrapolation warnings for all atoms are shown.
[in]resetewIf extrapolation warnings counter is reset every timestep.
[in]showewsumShow number of warnings every this many timesteps.
[in]maxewAbort simulation if more than this many warnings are encountered.
[in]cflengthLength unit conversion factor.
[in]cfenergyEnergy unit conversion factor.
[in]lammpsCutoffCutoff radius from LAMMPS (via pair_coeff).
[in]lammpsNtypesNumber of atom types in LAMMPS.
[in]myRankMPI process rank (passed on to structure index).

Definition at line 51 of file InterfaceLammps.cpp.

62{
63 this->emap = emap;
64 this->showew = showew;
65 this->resetew = resetew;
66 this->showewsum = showewsum;
67 this->maxew = maxew;
68 this->cflength = cflength;
69 this->cfenergy = cfenergy;
70 this->myRank = myRank;
71 log.writeToStdout = false;
72 string dir(directory);
73 char const separator = '/';
74 if (dir.back() != separator) dir += separator;
76 loadSettingsFile(dir + "input.nn");
77 setupGeneric(dir);
78 setupSymmetryFunctionScaling(dir + "scaling.data");
79 bool collectStatistics = false;
80 bool collectExtrapolationWarnings = false;
81 bool writeExtrapolationWarnings = false;
82 bool stopOnExtrapolationWarnings = false;
83 if (showew == true || showewsum > 0 || maxew >= 0)
84 {
85 collectExtrapolationWarnings = true;
86 }
87 setupSymmetryFunctionStatistics(collectStatistics,
88 collectExtrapolationWarnings,
90 stopOnExtrapolationWarnings);
92
93 log << "\n";
94 log << "*** SETUP: LAMMPS INTERFACE *************"
95 "**************************************\n";
96 log << "\n";
97
98 if (showew)
99 {
100 log << "Individual extrapolation warnings will be shown.\n";
101 }
102 else
103 {
104 log << "Individual extrapolation warnings will not be shown.\n";
105 }
106
107 if (showewsum != 0)
108 {
109 log << strpr("Extrapolation warning summary will be shown every %d"
110 " timesteps.\n", showewsum);
111 }
112 else
113 {
114 log << "Extrapolation warning summary will not be shown.\n";
115 }
116
117 if (maxew != 0)
118 {
119 log << strpr("The simulation will be stopped when %d extrapolation"
120 " warnings are exceeded.\n", maxew);
121 }
122 else
123 {
124 log << "No extrapolation warning limit set.\n";
125 }
126
127 if (resetew)
128 {
129 log << "Extrapolation warning counter is reset every time step.\n";
130 }
131 else
132 {
133 log << "Extrapolation warnings are accumulated over all time steps.\n";
134 }
135
136 log << "-----------------------------------------"
137 "--------------------------------------\n";
138 log << "CAUTION: If the LAMMPS unit system differs from the one used\n";
139 log << " during NN training, appropriate conversion factors\n";
140 log << " must be provided (see keywords cflength and cfenergy).\n";
141 log << "\n";
142 log << strpr("Length unit conversion factor: %24.16E\n", cflength);
143 log << strpr("Energy unit conversion factor: %24.16E\n", cfenergy);
144 double sfCutoff = getMaxCutoffRadius();
145 log << "\n";
146 log << "Checking consistency of cutoff radii (in LAMMPS units):\n";
147 log << strpr("LAMMPS Cutoff (via pair_coeff) : %11.3E\n", lammpsCutoff);
148 log << strpr("Maximum symmetry function cutoff: %11.3E\n", sfCutoff);
149 if (lammpsCutoff < sfCutoff)
150 {
151 throw runtime_error("ERROR: LAMMPS cutoff via pair_coeff keyword is"
152 " smaller than maximum symmetry function"
153 " cutoff.\n");
154 }
155 else if (fabs(sfCutoff - lammpsCutoff) / lammpsCutoff > TOLCUTOFF)
156 {
157 log << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!\n";
158 log << "WARNING: Potential length units mismatch!\n";
159 log << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!\n";
160 }
161 else
162 {
163 log << "Cutoff radii are consistent.\n";
164 }
165
166 log << "-----------------------------------------"
167 "--------------------------------------\n";
168 log << "Element mapping string from LAMMPS to n2p2: \""
169 + this->emap + "\"\n";
170 // Create default element mapping.
171 if (this->emap == "")
172 {
173 if (elementMap.size() != (size_t)lammpsNtypes)
174 {
175 throw runtime_error(strpr("ERROR: No element mapping given and "
176 "number of LAMMPS atom types (%d) and "
177 "NNP elements (%zu) does not match.\n",
178 lammpsNtypes, elementMap.size()));
179 }
180 log << "Element mapping string empty, creating default mapping.\n";
181 for (int i = 0; i < lammpsNtypes; ++i)
182 {
183 mapTypeToElement[i + 1] = i;
184 mapElementToType[i] = i + 1;
185 }
186 }
187 // Read element mapping from pair_style argument.
188 else
189 {
190 vector<string> emapSplit = split(reduce(trim(this->emap), " \t", ""),
191 ',');
192 if (elementMap.size() < emapSplit.size())
193 {
194 throw runtime_error(strpr("ERROR: Element mapping is inconsistent,"
195 " NNP elements: %zu,"
196 " emap elements: %zu.\n",
197 elementMap.size(),
198 emapSplit.size()));
199 }
200 for (string s : emapSplit)
201 {
202 vector<string> typeString = split(s, ':');
203 if (typeString.size() != 2)
204 {
205 throw runtime_error(strpr("ERROR: Invalid element mapping "
206 "string: \"%s\".\n", s.c_str()));
207 }
208 int t = stoi(typeString.at(0));
209 if (t > lammpsNtypes)
210 {
211 throw runtime_error(strpr("ERROR: LAMMPS type \"%d\" not "
212 "present, there are only %d types "
213 "defined.\n", t, lammpsNtypes));
214 }
215 size_t e = elementMap[typeString.at(1)];
216 mapTypeToElement[t] = e;
217 mapElementToType[e] = t;
218 }
219 }
220 log << "\n";
221 log << "CAUTION: Please ensure that this mapping between LAMMPS\n";
222 log << " atom types and NNP elements is consistent:\n";
223 log << "\n";
224 log << "---------------------------\n";
225 log << "LAMMPS type | NNP element\n";
226 log << "---------------------------\n";
227 for (int i = 1; i <= lammpsNtypes; ++i)
228 {
229 if (mapTypeToElement.find(i) != mapTypeToElement.end())
230 {
231 size_t e = mapTypeToElement.at(i);
232 log << strpr("%11d <-> %2s (%3zu)\n",
233 i,
234 elementMap[e].c_str(),
235 elementMap.atomicNumber(e));
236 ignoreType[i] = false;
237 }
238 else
239 {
240 log << strpr("%11d <-> --\n", i);
241 ignoreType[i] = true;
242
243 }
244 }
245 log << "---------------------------\n";
246 log << "\n";
247 log << "NNP setup for LAMMPS completed.\n";
248
249 log << "*****************************************"
250 "**************************************\n";
251
252 structure.setElementMap(elementMap);
253
254 initialized = true;
255}
#define TOLCUTOFF
Structure structure
Structure containing local atoms.
std::map< int, bool > ignoreType
True if atoms of this LAMMPS type will be ignored.
std::map< int, std::size_t > mapTypeToElement
Map from LAMMPS type to n2p2 element index.
void writeExtrapolationWarnings()
Write extrapolation warnings to log.
std::map< std::size_t, int > mapElementToType
Map from n2p2 element index to LAMMPS type.
double getMaxCutoffRadius() const
Get largest cutoff of symmetry functions.
std::string emap
Corresponds to LAMMPS map keyword.
ElementMap elementMap
Global element map, populated by setupElementMap().
Definition Mode.h:627
void initialize()
Write welcome message with version information.
Definition Mode.cpp:56
void setupGeneric(std::string const &nnpDir="", bool skipNormalize=false, bool initialHardness=false)
Combine multiple setup routines and provide a basic NNP setup.
Definition Mode.cpp:213
virtual void setupNeuralNetworkWeights(std::map< std::string, std::string > fileNameFormats=std::map< std::string, std::string >())
Set up neural network weights from files with given name format.
Definition Mode.cpp:1469
virtual void setupSymmetryFunctionScaling(std::string const &fileName="scaling.data")
Set up symmetry function scaling from file.
Definition Mode.cpp:736
Log log
Global log file.
Definition Mode.h:629
void setupSymmetryFunctionStatistics(bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
Set up symmetry function statistics collection.
Definition Mode.cpp:1127
void loadSettingsFile(std::string const &fileName="input.nn")
Open settings file and load all keywords into memory.
Definition Mode.cpp:162
string strpr(const char *format,...)
String version of printf function.
Definition utility.cpp:90
string trim(string const &line, string const &whitespace)
Remove leading and trailing whitespaces from string.
Definition utility.cpp:47
vector< string > split(string const &input, char delimiter)
Split string at each delimiter.
Definition utility.cpp:33
string reduce(string const &line, string const &whitespace, string const &fill)
Replace multiple whitespaces with fill.
Definition utility.cpp:60

References cfenergy, cflength, nnp::Mode::elementMap, emap, getMaxCutoffRadius(), ignoreType, nnp::Mode::initialize(), initialized, nnp::Mode::loadSettingsFile(), nnp::Mode::log, mapElementToType, mapTypeToElement, maxew, myRank, nnp::reduce(), resetew, nnp::Mode::setupGeneric(), nnp::Mode::setupNeuralNetworkWeights(), nnp::Mode::setupSymmetryFunctionScaling(), nnp::Mode::setupSymmetryFunctionStatistics(), showew, showewsum, nnp::split(), nnp::strpr(), structure, TOLCUTOFF, nnp::trim(), and writeExtrapolationWarnings().

Here is the call graph for this function:

◆ setGlobalStructureStatus()

void InterfaceLammps::setGlobalStructureStatus ( bool const status)

Specify whether n2p2 knows about global structure or only local structure.

Parameters
[in]statustrue if n2p2 has global structure.

Definition at line 257 of file InterfaceLammps.cpp.

258{
259 hasGlobalStructure = status;
260}

References hasGlobalStructure.

◆ getGlobalStructureStatus()

bool InterfaceLammps::getGlobalStructureStatus ( )

Check if n2p2 knows about global structure.

Definition at line 262 of file InterfaceLammps.cpp.

263{
264 return hasGlobalStructure;
265}

References hasGlobalStructure.

◆ setLocalAtoms()

void InterfaceLammps::setLocalAtoms ( int numAtomsLocal,
int const *const atomType )

(Re)set structure to contain only local LAMMPS atoms.

Parameters
[in]numAtomsLocalNumber of local atoms.
[in]atomTypeLAMMPS atom type.

Definition at line 267 of file InterfaceLammps.cpp.

269{
270 for (size_t i = 0; i < numElements; ++i)
271 {
272 structure.numAtomsPerElement[i] = 0;
273 }
274 structure.index = myRank;
275 structure.numAtoms = 0;
276 structure.hasNeighborList = false;
277 structure.hasSymmetryFunctions = false;
278 structure.hasSymmetryFunctionDerivatives = false;
279 structure.energy = 0.0;
280 structure.atoms.clear();
281 indexMap.clear();
282 structure.atoms.reserve(numAtomsLocal);
283 indexMap.resize(numAtomsLocal, numeric_limits<size_t>::max());
284 for (int i = 0; i < numAtomsLocal; i++)
285 {
286 if (ignoreType[atomType[i]]) continue;
287 indexMap.at(i) = structure.numAtoms;
288 structure.numAtoms++;
289 structure.atoms.push_back(Atom());
290 Atom& a = structure.atoms.back();
291 a.index = i;
293 a.element = mapTypeToElement[atomType[i]];
294 a.numNeighbors = 0;
295 a.hasSymmetryFunctions = false;
297 a.neighbors.clear();
298 a.numNeighborsPerElement.clear();
300 structure.numAtomsPerElement[a.element]++;
301 }
302
303 return;
304}
std::vector< size_t > indexMap
Map from LAMMPS index to n2p2 atom index.
std::size_t numElements
Definition Mode.h:667
std::vector< Neighbor > neighbors
Neighbor array (maximum number defined in macros.h.
Definition Atom.h:170
bool hasSymmetryFunctionDerivatives
If symmetry function derivatives are saved for this atom.
Definition Atom.h:97
std::size_t index
Index number of this atom.
Definition Atom.h:101
std::size_t indexStructure
Index number of structure this atom belongs to.
Definition Atom.h:103
std::size_t element
Element index of this atom.
Definition Atom.h:107
bool hasSymmetryFunctions
If symmetry function values are saved for this atom.
Definition Atom.h:95
std::vector< std::size_t > numNeighborsPerElement
Number of neighbors per element.
Definition Atom.h:138
std::size_t numNeighbors
Total number of neighbors.
Definition Atom.h:109

References nnp::Atom::element, nnp::Atom::hasSymmetryFunctionDerivatives, nnp::Atom::hasSymmetryFunctions, ignoreType, nnp::Atom::index, indexMap, nnp::Atom::indexStructure, mapTypeToElement, myRank, nnp::Atom::neighbors, nnp::Mode::numElements, nnp::Atom::numNeighbors, nnp::Atom::numNeighborsPerElement, and structure.

◆ setLocalAtomPositions()

void InterfaceLammps::setLocalAtomPositions ( double const *const *const atomPos)

Set absolute atom positions from LAMMPS (nnp/develop only).

Parameters
[in]atomPosAtom coordinate array in LAMMPS units.

Definition at line 306 of file InterfaceLammps.cpp.

307{
308 for (size_t i = 0; i < structure.numAtoms; ++i)
309 {
310 Atom& a = structure.atoms.at(i);
311 a.r[0] = atomPos[i][0] * cflength;
312 a.r[1] = atomPos[i][1] * cflength;
313 a.r[2] = atomPos[i][2] * cflength;
314 if (normalize)
315 {
316 a.r[0] *= convLength;
317 a.r[1] *= convLength;
318 a.r[2] *= convLength;
319 }
320 }
321
322 return;
323}
bool normalize
Definition Mode.h:665
double convLength
Definition Mode.h:674
Vec3D r
Cartesian coordinates.
Definition Atom.h:125

References cflength, nnp::Mode::convLength, nnp::Mode::normalize, nnp::Atom::r, and structure.

◆ setLocalTags() [1/2]

void InterfaceLammps::setLocalTags ( int const *const atomTag)

Set atom tags (int version, -DLAMMPS_SMALLBIG).

Parameters
[in]atomTagLAMMPS atom tag.

Definition at line 325 of file InterfaceLammps.cpp.

326{
327 for (size_t i = 0; i < structure.atoms.size(); i++)
328 {
329 // Implicit conversion from int to int64_t!
330 structure.atoms.at(i).tag = atomTag[i];
331 }
332
333 return;
334}

References structure.

◆ setLocalTags() [2/2]

void InterfaceLammps::setLocalTags ( int64_t const *const atomTag)

Set atom tags (int64_t version, -DLAMMPS_BIGBIG).

Parameters
[in]atomTagLAMMPS atom tag.

Definition at line 336 of file InterfaceLammps.cpp.

337{
338 for (size_t i = 0; i < structure.atoms.size(); i++)
339 {
340 structure.atoms.at(i).tag = atomTag[i];
341 }
342
343 return;
344}

References structure.

◆ setBoxVectors()

void InterfaceLammps::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).

Parameters
[in]boxloArray containing coordinates of origin xlo, ylo, zlo.
[in]boxhiArray containing coordinates xhi, yhi, zhi.
[in]xyTilt factor for box vector b.
[in]xzFirst tilt factor for box vector c.
[in]yzSecond tilt factor for box vector c.

Definition at line 346 of file InterfaceLammps.cpp.

351{
352 structure.isPeriodic = true;
353
354 // Box vector a
355 structure.box[0][0] = boxhi[0] - boxlo[0];
356 structure.box[0][1] = 0;
357 structure.box[0][2] = 0;
358
359 // Box vector b
360 structure.box[1][0] = xy;
361 structure.box[1][1] = boxhi[1] - boxlo[1];
362 structure.box[1][2] = 0;
363
364 // Box vector c
365 structure.box[2][0] = xz;
366 structure.box[2][1] = yz;
367 structure.box[2][2] = boxhi[2] - boxlo[2];
368
369 // LAMMPS may set triclinic = 1 even if the following condition is not
370 // satisfied.
371 if (structure.box[0][1] > numeric_limits<double>::min() ||
372 structure.box[0][2] > numeric_limits<double>::min() ||
373 structure.box[1][0] > numeric_limits<double>::min() ||
374 structure.box[1][2] > numeric_limits<double>::min() ||
375 structure.box[2][0] > numeric_limits<double>::min() ||
376 structure.box[2][1] > numeric_limits<double>::min())
377 {
378 structure.isTriclinic = true;
379 }
380
381 for(size_t i = 0; i < 3; ++i)
382 {
383 structure.box[i] *= cflength;
384 if (normalize) structure.box[i] *= convLength;
385 }
386
387 structure.calculateInverseBox();
388 structure.calculateVolume();
389 //cout << "Box vectors: \n";
390 //for(size_t i = 0; i < 3; ++i)
391 //{
392 // for(size_t j = 0; j < 3; ++j)
393 // {
394 // cout << structure.box[i][j] / convLength << " ";
395 // }
396 // cout << endl;
397 //}
398
399}

References cflength, nnp::Mode::convLength, nnp::Mode::normalize, and structure.

◆ allocateNeighborlists()

void InterfaceLammps::allocateNeighborlists ( int const *const numneigh)

Allocate neighbor lists.

Parameters
[in]numneighArray containing number of neighbors for each local atom.

Definition at line 401 of file InterfaceLammps.cpp.

402{
403 for(size_t i = 0; i < structure.numAtoms; ++i)
404 {
405 auto& a = structure.atoms.at(i);
406 a.neighbors.reserve(numneigh[i]);
407 }
408}

References structure.

◆ addNeighbor()

void InterfaceLammps::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).

Parameters
[in]iLocal atom index.
[in]jNeighbor atom index.
[in]tagNeighbor atom tag.
[in]typeNeighbor atom type.
[in]dxNeighbor atom distance in x direction.
[in]dyNeighbor atom distance in y direction.
[in]dzNeighbor atom distance in z direction.
[in]d2Square of neighbor atom distance.

If -DLAMMPS_SMALLBIG implicit conversion is applied for tag.

Definition at line 410 of file InterfaceLammps.cpp.

418{
419 if (ignoreType[type] ||
420 indexMap.at(i) == numeric_limits<size_t>::max()) return;
421 Atom& a = structure.atoms[indexMap.at(i)];
422 a.numNeighbors++;
423 a.neighbors.push_back(Atom::Neighbor());
425 Atom::Neighbor& n = a.neighbors.back();
426
427 n.index = j;
428 n.tag = tag;
429 n.element = mapTypeToElement[type];
430 n.dr[0] = dx * cflength;
431 n.dr[1] = dy * cflength;
432 n.dr[2] = dz * cflength;
433 n.d = sqrt(d2) * cflength;
434 if (normalize)
435 {
436 n.dr *= convLength;
437 n.d *= convLength;
438 }
439
440 return;
441}
std::size_t index
Index of neighbor atom.
Definition Atom.h:38
std::size_t element
Element index of neighbor atom.
Definition Atom.h:42
double d
Distance to neighbor atom.
Definition Atom.h:44
Vec3D dr
Distance vector to neighbor atom.
Definition Atom.h:46
int64_t tag
Tag of neighbor atom.
Definition Atom.h:40

References cflength, nnp::Mode::convLength, nnp::Atom::Neighbor::d, nnp::Atom::Neighbor::dr, nnp::Atom::Neighbor::element, ignoreType, nnp::Atom::Neighbor::index, indexMap, mapTypeToElement, nnp::Atom::neighbors, nnp::Mode::normalize, nnp::Atom::numNeighbors, nnp::Atom::numNeighborsPerElement, structure, and nnp::Atom::Neighbor::tag.

◆ finalizeNeighborList()

void InterfaceLammps::finalizeNeighborList ( )

Sorts neighbor list and creates cutoff map if necessary.

If structure is periodic, this function needs to be called after setBoxVectors!

Definition at line 443 of file InterfaceLammps.cpp.

444{
446 {
447 for (auto& a : structure.atoms)
448 {
449 a.hasNeighborList = true;
450 }
451 // Ewald summation cut-off depends on box vectors.
452 structure.calculateMaxCutoffRadiusOverall(
454 screeningFunction.getOuter(),
456 structure.sortNeighborList();
457 structure.setupNeighborCutoffMap(cutoffs);
458 }
459
460 return;
461}
std::vector< std::vector< double > > cutoffs
Matrix storing all symmetry function cut-offs for all elements.
Definition Mode.h:689
NNPType nnpType
Definition Mode.h:664
double maxCutoffRadius
Definition Mode.h:670
@ HDNNP_4G
NNP with electrostatics and non-local charge transfer (4G-HDNNP).
Definition Mode.h:95
ScreeningFunction screeningFunction
Definition Mode.h:682
EwaldSetup ewaldSetup
Definition Mode.h:677

References nnp::Mode::cutoffs, nnp::Mode::ewaldSetup, nnp::Mode::HDNNP_4G, nnp::Mode::maxCutoffRadius, nnp::Mode::nnpType, nnp::Mode::screeningFunction, and structure.

◆ process()

void InterfaceLammps::process ( )

Calculate symmetry functions, atomic neural networks and sum of local energy contributions.

Definition at line 463 of file InterfaceLammps.cpp.

464{
465#ifdef N2P2_NO_SF_GROUPS
467#else
469#endif
471 {
474 if (normalize)
475 {
476 structure.energy = physicalEnergy(structure, false);
477 }
479 }
480 else if (nnpType == NNPType::HDNNP_4G)
481 {
482 if (!isElecDone)
483 {
485 isElecDone = true;
486 } else
487 {
489 isElecDone = false;
491 if (normalize)
492 {
493 structure.energy = physicalEnergy(structure, false);
494 }
496 }
497 }
498
499 return;
500}
double physicalEnergy(Structure const &structure, bool ref=true) const
Undo normalization for a given energy of structure.
Definition Mode.cpp:2138
void addEnergyOffset(Structure &structure, bool ref=true)
Add atomic energy offsets to reference energy.
Definition Mode.cpp:2034
@ HDNNP_2G
Short range NNP (2G-HDNNP).
Definition Mode.h:93
void calculateAtomicNeuralNetworks(Structure &structure, bool const derivatives, std::string id="")
Calculate a single atomic neural network for a given atom and nn type.
Definition Mode.cpp:1666
void calculateEnergy(Structure &structure) const
Calculate potential energy for a given structure.
Definition Mode.cpp:1831
void calculateSymmetryFunctionGroups(Structure &structure, bool const derivatives)
Calculate all symmetry function groups for all atoms in given structure.
Definition Mode.cpp:1585
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
Definition Mode.cpp:1504

References nnp::Mode::addEnergyOffset(), nnp::Mode::calculateAtomicNeuralNetworks(), nnp::Mode::calculateEnergy(), nnp::Mode::calculateSymmetryFunctionGroups(), nnp::Mode::calculateSymmetryFunctions(), nnp::Mode::HDNNP_2G, nnp::Mode::HDNNP_4G, isElecDone, nnp::Mode::nnpType, nnp::Mode::normalize, nnp::Mode::physicalEnergy(), and structure.

Here is the call graph for this function:

◆ processDevelop()

void InterfaceLammps::processDevelop ( )

Calculate symmetry functions, atomic neural networks and sum of local energy contributions (development version for "hdnnp/develop" pair style).

Definition at line 502 of file InterfaceLammps.cpp.

503{
504#ifdef N2P2_NO_SF_GROUPS
506#else
508#endif
511 {
514 ewaldSetup.logEwaldCutoffs(log, convLength * cflength);
515 }
517 if (nnpType == NNPType::HDNNP_4G ||
519 if (normalize)
520 {
521 structure.energy = physicalEnergy(structure, false);
522 }
524
525 return;
526}
@ HDNNP_Q
Short range NNP with charge NN, no electrostatics/Qeq (M.
Definition Mode.h:103
void calculateCharge(Structure &structure) const
Calculate total charge for a given structure.
Definition Mode.cpp:1858
void chargeEquilibration(Structure &structure, bool const derivativesElec)
Perform global charge equilibration method.
Definition Mode.cpp:1779

References nnp::Mode::addEnergyOffset(), nnp::Mode::calculateAtomicNeuralNetworks(), nnp::Mode::calculateCharge(), nnp::Mode::calculateEnergy(), nnp::Mode::calculateSymmetryFunctionGroups(), nnp::Mode::calculateSymmetryFunctions(), cflength, nnp::Mode::chargeEquilibration(), nnp::Mode::convLength, nnp::Mode::ewaldSetup, nnp::Mode::HDNNP_4G, nnp::Mode::HDNNP_Q, nnp::Mode::log, nnp::Mode::nnpType, nnp::Mode::normalize, nnp::Mode::physicalEnergy(), and structure.

Here is the call graph for this function:

◆ getEnergy()

double InterfaceLammps::getEnergy ( ) const

Return sum of local energy contributions.

Returns
Sum of local energy contributions.

Definition at line 550 of file InterfaceLammps.cpp.

551{
552 return structure.energy / cfenergy;
553}

References cfenergy, and structure.

◆ getAtomicEnergy()

double InterfaceLammps::getAtomicEnergy ( int index) const

Return energy contribution of one atom.

Parameters
[in]indexAtom index.
Returns
energy contribution of atom with given index.
Attention
These atomic contributions are not physical!

Definition at line 555 of file InterfaceLammps.cpp.

556{
557 Atom const& a = structure.atoms.at(index);
558 Element const& e = elements.at(a.element);
559
560 if (normalize)
561 {
562 return (physical("energy", a.energy)
563 + meanEnergy
565 }
566 else
567 {
568 return (a.energy + e.getAtomicEnergyOffset()) / cfenergy;
569 }
570}
double getAtomicEnergyOffset() const
Get atomicEnergyOffset.
Definition Element.h:315
double meanEnergy
Definition Mode.h:672
double physical(std::string const &property, double value) const
Undo normalization for a given property.
Definition Mode.cpp:2126
std::vector< Element > elements
Definition Mode.h:683
double energy
Atomic energy determined by neural network.
Definition Atom.h:115

References cfenergy, nnp::Atom::element, nnp::Mode::elements, nnp::Atom::energy, nnp::Element::getAtomicEnergyOffset(), nnp::Mode::meanEnergy, nnp::Mode::normalize, nnp::Mode::physical(), and structure.

Here is the call graph for this function:

◆ addElectrostaticEnergy()

void InterfaceLammps::addElectrostaticEnergy ( double energy)

Adds electrostatic energy contribution to the total structure energy.

Parameters
[in]electrostaticenergy (calculated in LAMMPS).

Definition at line 674 of file InterfaceLammps.cpp.

675{
676 structure.energyElec = eElec;
677}

References structure.

◆ getForces()

void InterfaceLammps::getForces ( double *const *const & atomF) const

Calculate forces and add to LAMMPS atomic force arrays.

Parameters
[in,out]atomFLAMMPS force array for local and ghost atoms.

Definition at line 679 of file InterfaceLammps.cpp.

679 {
680 double const cfforce = cflength / cfenergy;
681 double convForce = 1.0;
682 if (normalize) {
683 convForce = convLength / convEnergy;
684 }
685
686 // Loop over all local atoms. Neural network and Symmetry function
687 // derivatives are saved in the dEdG arrays of atoms and dGdr arrays of
688 // atoms and their neighbors. These are now summed up to the force
689 // contributions of local and ghost atoms.
690 Atom const *a = NULL;
691
692 for (size_t i = 0; i < structure.atoms.size(); ++i) {
693 // Set pointer to atom.
694 a = &(structure.atoms.at(i));
695
696#ifndef NNP_FULL_SFD_MEMORY
697 vector <vector<size_t>> const &tableFull
698 = elements.at(a->element).getSymmetryFunctionTable();
699#endif
700 // Loop over all neighbor atoms. Some are local, some are ghost atoms.
701 for (vector<Atom::Neighbor>::const_iterator n = a->neighbors.begin();
702 n != a->neighbors.end(); ++n) {
703 // Temporarily save the neighbor index. Note: this is the index for
704 // the LAMMPS force array.
705 size_t const in = n->index;
706 // Now loop over all symmetry functions and add force contributions
707 // (local + ghost atoms).
708#ifndef NNP_FULL_SFD_MEMORY
709 vector <size_t> const &table = tableFull.at(n->element);
710 for (size_t s = 0; s < n->dGdr.size(); ++s)
711 {
712 double const dEdG = a->dEdG[table.at(s)] * cfforce * convForce;
713#else
714 for (size_t s = 0; s < a->numSymmetryFunctions; ++s)
715 {
716 double const dEdG = a->dEdG[s] * cfforce * convForce;
717#endif
718 double const *const dGdr = n->dGdr[s].r;
719 atomF[in][0] -= dEdG * dGdr[0];
720 atomF[in][1] -= dEdG * dGdr[1];
721 atomF[in][2] -= dEdG * dGdr[2];
722 }
723 }
724 // Temporarily save the atom index. Note: this is the index for
725 // the LAMMPS force array.
726 size_t const ia = a->index;
727 // Loop over all symmetry functions and add force contributions (local
728 // atoms).
729 for (size_t s = 0; s < a->numSymmetryFunctions; ++s)
730 {
731 double const dEdG = a->dEdG[s] * cfforce * convForce;
732 double const *const dGdr = a->dGdr[s].r;
733 atomF[ia][0] -= dEdG * dGdr[0];
734 atomF[ia][1] -= dEdG * dGdr[1];
735 atomF[ia][2] -= dEdG * dGdr[2];
736 }
737 }
738
739 return;
740}
double convEnergy
Definition Mode.h:673
std::size_t numSymmetryFunctions
Number of symmetry functions used to describe the atom environment.
Definition Atom.h:113
std::vector< double > dEdG
Derivative of atomic energy with respect to symmetry functions.
Definition Atom.h:149
std::vector< Vec3D > dGdr
Derivative of symmetry functions with respect to this atom's coordinates.
Definition Atom.h:161

References cfenergy, cflength, nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Atom::dEdG, nnp::Atom::dGdr, nnp::Atom::element, nnp::Mode::elements, nnp::Atom::index, nnp::Atom::neighbors, nnp::Mode::normalize, nnp::Atom::numSymmetryFunctions, and structure.

◆ getForcesDevelop()

void InterfaceLammps::getForcesDevelop ( double *const *const & atomF) const

Calculate forces and add to LAMMPS atomic force arrays (development version for "hdnnp/develop" pair style).

Parameters
[in,out]atomFLAMMPS force array for local and ghost atoms.

Definition at line 742 of file InterfaceLammps.cpp.

743{
744 double const cfforce = cflength / cfenergy;
745 double convForce = 1.0;
746 if (normalize)
747 {
748 convForce = convLength / convEnergy;
749 }
750
751 // Loop over all local atoms. Neural network and Symmetry function
752 // derivatives are saved in the dEdG arrays of atoms and dGdr arrays of
753 // atoms and their neighbors. These are now summed up to the force
754 // contributions of local and ghost atoms.
755 for (auto const& a : structure.atoms)
756 {
757 size_t const ia = a.index;
758 Vec3D selfForce = a.calculateSelfForceShort();
759 selfForce *= cfforce * convForce;
760 // TODO: ia is not the right index when some atom types are excluded / ignored
761 // (see use of indexmap)
762 add3DVecToArray(atomF[ia], selfForce);
763
764#ifndef N2P2_FULL_SFD_MEMORY
765 vector<vector<size_t> > const& tableFull
766 = elements.at(a.element).getSymmetryFunctionTable();
767#endif
768 // Loop over all neighbor atoms. Some are local, some are ghost atoms.
769
770 //for (auto const& n : a.neighbors)
771 size_t const numNeighbors = a.getStoredMinNumNeighbors(maxCutoffRadius);
772#ifdef _OPENMP
773 #pragma omp parallel for
774#endif
775 for (size_t k = 0; k < numNeighbors; ++k)
776 {
777 Atom::Neighbor const& n = a.neighbors[k];
778 // Temporarily save the neighbor index. Note: this is the index for
779 // the LAMMPS force array.
780 size_t const in = n.index;
781
782#ifndef N2P2_FULL_SFD_MEMORY
783 Vec3D pairForce = a.calculatePairForceShort(n, &tableFull);
784#else
785 Vec3D pairForce = a.calculatePairForceShort(n);
786#endif
787 pairForce *= cfforce * convForce;
788 add3DVecToArray(atomF[in], pairForce);
789 }
790 }
791
792 // Comment: Will not work with multiple MPI tasks but this routine will
793 // probably be obsolete when Emir's solution is finished.
795 {
796 Structure const& s = structure;
797 VectorXd lambdaTotal = s.calculateForceLambdaTotal();
798
799#ifdef _OPENMP
800 #pragma omp parallel for
801#endif
802 // OpenMP 4.0 doesn't support range based loops
803 for (size_t i = 0; i < s.numAtoms; ++i)
804 {
805 auto const& ai = s.atoms[i];
806 add3DVecToArray(atomF[i], -ai.pEelecpr * cfforce * convForce);
807
808 for (auto const& aj : s.atoms)
809 {
810 size_t const j = aj.index;
811
812#ifndef N2P2_FULL_SFD_MEMORY
813 vector<vector<size_t> > const& tableFull
814 = elements.at(aj.element).getSymmetryFunctionTable();
815 Vec3D dChidr = aj.calculateDChidr(ai.index,
817 &tableFull);
818#else
819 Vec3D dChidr = aj.calculateDChidr(ai.index,
821#endif
822
823 Vec3D remainingForce = -lambdaTotal(j) * (ai.dAdrQ[j] + dChidr);
824 add3DVecToArray(atomF[i], remainingForce * cfforce * convForce);
825
826 }
827 }
828 }
829 return;
830}
void add3DVecToArray(double *const &arr, Vec3D const &v) const
Add a Vec3D vector to a 3D array in place.
Eigen::VectorXd const calculateForceLambdaTotal() const
Calculate lambda_total vector which is needed for the total force calculation in 4G NN.
std::size_t numAtoms
Total number of atoms present in this structure.
Definition Structure.h:75
std::vector< Atom > atoms
Vector of all atoms in this structure.
Definition Structure.h:122

References add3DVecToArray(), nnp::Structure::atoms, nnp::Structure::calculateForceLambdaTotal(), nnp::Atom::calculatePairForceShort(), nnp::Atom::calculateSelfForceShort(), cfenergy, cflength, nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Atom::element, nnp::Mode::elements, nnp::Atom::getStoredMinNumNeighbors(), nnp::Mode::HDNNP_4G, nnp::Atom::index, nnp::Atom::Neighbor::index, nnp::Mode::maxCutoffRadius, nnp::Atom::neighbors, nnp::Mode::nnpType, nnp::Mode::normalize, nnp::Structure::numAtoms, and structure.

Here is the call graph for this function:

◆ getForcesChi()

void InterfaceLammps::getForcesChi ( double const *const & lambda,
double *const *const & atomF ) const

Calculate chi-term for forces and add to LAMMPS atomic force arrays.

Parameters
[in,out]atomFLAMMPS force array for local and ghost atoms.

Definition at line 832 of file InterfaceLammps.cpp.

833 {
834 double const cfforce = cflength / cfenergy;
835 double convForce = 1.0;
836 if (normalize) {
837 convForce = convLength / convEnergy;
838 }
839
840 // Loop over all local atoms. Neural network and Symmetry function
841 // derivatives are saved in the dEdG arrays of atoms and dGdr arrays of
842 // atoms and their neighbors. These are now summed up to the force
843 // contributions of local and ghost atoms.
844 Atom const *a = NULL;
845
846 for (size_t i = 0; i < structure.atoms.size(); ++i) {
847 // Set pointer to atom.
848 a = &(structure.atoms.at(i));
849
850 // Temporarily save the atom index. Note: this is the index for
851 // the LAMMPS force array.
852 size_t const ia = a->index;
853 // Also save tag - 1 which is the correct position in lambda array.
854 size_t const ta = a->tag - 1;
855
856#ifndef NNP_FULL_SFD_MEMORY
857 vector <vector<size_t>> const &tableFull
858 = elements.at(a->element).getSymmetryFunctionTable();
859#endif
860 // Loop over all neighbor atoms. Some are local, some are ghost atoms.
861 for (vector<Atom::Neighbor>::const_iterator n = a->neighbors.begin();
862 n != a->neighbors.end(); ++n) {
863 // Temporarily save the neighbor index. Note: this is the index for
864 // the LAMMPS force array.
865 size_t const in = n->index;
866 // Also save tag - 1 which is the correct position in lambda array.
867 size_t const tn = n->tag - 1;
868 //std::cout << "Chi : " << a->chi << '\t' << "nei :" << n->index << '\n';
869 // Now loop over all symmetry functions and add force contributions
870 // (local + ghost atoms).
871#ifndef NNP_FULL_SFD_MEMORY
872 vector <size_t> const &table = tableFull.at(n->element);
873 for (size_t s = 0; s < n->dGdr.size(); ++s)
874 {
875 double const dChidG = a->dChidG[table.at(s)]
876 * cfforce * convForce;
877#else
878 for (size_t s = 0; s < a->numSymmetryFunctions; ++s)
879 {
880 double const dChidG = a->dChidG[s] * cfforce * convForce;
881#endif
882 double const *const dGdr = n->dGdr[s].r;
883 atomF[in][0] -= lambda[ta] * dChidG * dGdr[0];
884 atomF[in][1] -= lambda[ta] * dChidG * dGdr[1];
885 atomF[in][2] -= lambda[ta] * dChidG * dGdr[2];
886 }
887 }
888 // Loop over all symmetry functions and add force contributions (local
889 // atoms).
890 for (size_t s = 0; s < a->numSymmetryFunctions; ++s)
891 {
892 double const dChidG = a->dChidG[s] * cfforce * convForce;
893 double const *const dGdr = a->dGdr[s].r;
894 atomF[ia][0] -= lambda[ta] * dChidG * dGdr[0];
895 atomF[ia][1] -= lambda[ta] * dChidG * dGdr[1];
896 atomF[ia][2] -= lambda[ta] * dChidG * dGdr[2];
897 }
898 }
899
900 return;
901}
int64_t tag
Tag number of this atom.
Definition Atom.h:105
std::vector< double > dChidG
Derivative of electronegativity with respect to symmetry functions.
Definition Atom.h:153

References cfenergy, cflength, nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Atom::dChidG, nnp::Atom::dGdr, nnp::Atom::element, nnp::Mode::elements, nnp::Atom::index, nnp::Atom::neighbors, nnp::Mode::normalize, nnp::Atom::numSymmetryFunctions, structure, and nnp::Atom::tag.

◆ getCharges()

void InterfaceLammps::getCharges ( double *const & atomQ) const

Transfer charges (in units of e) to LAMMPS atomic charge vector.

Call after getAtomicEnergy().

Parameters
[in,out]atomQLAMMPS charge vector.

Definition at line 903 of file InterfaceLammps.cpp.

904{
905 if (nnpType != NNPType::HDNNP_4G) return;
906 if (!atomQ) return;
907
908 Structure const& s = structure;
909#ifdef _OPENMP
910 #pragma omp parallel for
911#endif
912 for (size_t i = 0; i < s.numAtoms; ++i)
913 {
914 atomQ[i] = s.atoms[i].charge;
915 }
916
917 return;
918}

References nnp::Structure::atoms, nnp::Mode::HDNNP_4G, nnp::Mode::nnpType, nnp::Structure::numAtoms, and structure.

◆ isInitialized()

bool nnp::InterfaceLammps::isInitialized ( ) const
inline

Check if this interface is correctly initialized.

Returns
True if initialized, False otherwise.

Definition at line 326 of file InterfaceLammps.h.

327{
328 return initialized;
329}

References initialized.

◆ getMaxCutoffRadius()

double InterfaceLammps::getMaxCutoffRadius ( ) const

Get largest cutoff of symmetry functions.

Returns
Largest cutoff of all symmetry functions.

Definition at line 528 of file InterfaceLammps.cpp.

529{
531 else return maxCutoffRadius / cflength;
532}

References cflength, nnp::Mode::convLength, nnp::Mode::maxCutoffRadius, and nnp::Mode::normalize.

Referenced by getMaxCutoffRadiusOverall(), and initialize().

Here is the caller graph for this function:

◆ getEwaldPrec()

double InterfaceLammps::getEwaldPrec ( ) const

Get Ewald precision parameter.

Returns
Ewald precision parameter.

Definition at line 664 of file InterfaceLammps.cpp.

665{
667}
double getEwaldPrecision() const
Getter for Mode::ewaldSetup.precision.
Definition Mode.h:740

References nnp::Mode::getEwaldPrecision().

Here is the call graph for this function:

◆ getMaxCutoffRadiusOverall()

double InterfaceLammps::getMaxCutoffRadiusOverall ( )

Get largest cutoff including structure specific cutoff and screening cutoff.

Returns
Largest cutoff of all symmetry functions and structure specific cutoff and screening cutoff.

Definition at line 534 of file InterfaceLammps.cpp.

535{
536 double cutoff = 0;
538 {
539 structure.calculateMaxCutoffRadiusOverall(
541 screeningFunction.getOuter(),
543 cutoff = structure.maxCutoffRadiusOverall / cflength;
544 if (normalize) cutoff /= convLength;
545 }
546 else cutoff = getMaxCutoffRadius();
547 return cutoff;
548}

References cflength, nnp::Mode::convLength, nnp::Mode::ewaldSetup, getMaxCutoffRadius(), nnp::Mode::HDNNP_4G, nnp::Mode::maxCutoffRadius, nnp::Mode::nnpType, nnp::Mode::normalize, nnp::Mode::screeningFunction, and structure.

Here is the call graph for this function:

◆ getEWBufferSize()

long InterfaceLammps::getEWBufferSize ( ) const

Calculate buffer size for extrapolation warning communication.

Returns
Buffer size.

Definition at line 920 of file InterfaceLammps.cpp.

921{
922 long bs = 0;
923#ifndef N2P2_NO_MPI
924 int ss = 0; // size_t size.
925 int ds = 0; // double size.
926 int cs = 0; // char size.
927 MPI_Pack_size(1, MPI_SIZE_T, MPI_COMM_WORLD, &ss);
928 MPI_Pack_size(1, MPI_DOUBLE, MPI_COMM_WORLD, &ds);
929 MPI_Pack_size(1, MPI_CHAR , MPI_COMM_WORLD, &cs);
930
931 for (vector<Element>::const_iterator it = elements.begin();
932 it != elements.end(); ++it)
933 {
934 map<size_t, SymFncStatistics::Container> const& m
935 = it->statistics.data;
936 bs += ss; // n.
937 for (map<size_t, SymFncStatistics::Container>::const_iterator
938 it2 = m.begin(); it2 != m.end(); ++it2)
939 {
940 bs += ss; // index (it2->first).
941 bs += ss; // countEW (it2->second.countEW).
942 bs += ss; // type (it2->second.type).
943 bs += ds; // Gmin (it2->second.Gmin).
944 bs += ds; // Gmax (it2->second.Gmax).
945 bs += ss; // element.length() (it2->second.element.length()).
946 bs += (it2->second.element.length() + 1) * cs; // element.
947 size_t countEW = it2->second.countEW;
948 bs += countEW * ss; // indexStructureEW.
949 bs += countEW * ss; // indexAtomEW.
950 bs += countEW * ds; // valueEW.
951 }
952 }
953#endif
954 return bs;
955}
#define MPI_SIZE_T
Definition mpi-extra.h:22

References nnp::Mode::elements, and MPI_SIZE_T.

◆ fillEWBuffer()

void InterfaceLammps::fillEWBuffer ( char *const & buf,
int bs ) const

Fill provided buffer with extrapolation warning entries.

Parameters
[in,out]bufCommunication buffer to fill.
[in]bsBuffer size.

Definition at line 957 of file InterfaceLammps.cpp.

958{
959#ifndef N2P2_NO_MPI
960 int p = 0;
961 for (vector<Element>::const_iterator it = elements.begin();
962 it != elements.end(); ++it)
963 {
964 map<size_t, SymFncStatistics::Container> const& m =
965 it->statistics.data;
966 size_t n = m.size();
967 MPI_Pack((void *) &(n), 1, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
968 for (map<size_t, SymFncStatistics::Container>::const_iterator
969 it2 = m.begin(); it2 != m.end(); ++it2)
970 {
971 MPI_Pack((void *) &(it2->first ), 1, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
972 size_t countEW = it2->second.countEW;
973 MPI_Pack((void *) &(countEW ), 1, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
974 MPI_Pack((void *) &(it2->second.type ), 1, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
975 MPI_Pack((void *) &(it2->second.Gmin ), 1, MPI_DOUBLE, buf, bs, &p, MPI_COMM_WORLD);
976 MPI_Pack((void *) &(it2->second.Gmax ), 1, MPI_DOUBLE, buf, bs, &p, MPI_COMM_WORLD);
977 // it2->element
978 size_t ts = it2->second.element.length() + 1;
979 MPI_Pack((void *) &ts , 1, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
980 MPI_Pack((void *) it2->second.element.c_str() , ts, MPI_CHAR , buf, bs, &p, MPI_COMM_WORLD);
981 MPI_Pack((void *) &(it2->second.indexStructureEW.front()), countEW, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
982 MPI_Pack((void *) &(it2->second.indexAtomEW.front() ), countEW, MPI_SIZE_T, buf, bs, &p, MPI_COMM_WORLD);
983 MPI_Pack((void *) &(it2->second.valueEW.front() ), countEW, MPI_DOUBLE, buf, bs, &p, MPI_COMM_WORLD);
984 }
985 }
986#endif
987 return;
988}
size_t p

References nnp::Mode::elements, MPI_SIZE_T, and p.

◆ extractEWBuffer()

void InterfaceLammps::extractEWBuffer ( char const *const & buf,
int bs )

Extract given buffer to symmetry function statistics class.

Parameters
[in]bufBuffer with extrapolation warnings data.
[in]bsBuffer size.

Definition at line 990 of file InterfaceLammps.cpp.

991{
992#ifndef N2P2_NO_MPI
993 int p = 0;
994 for (vector<Element>::iterator it = elements.begin();
995 it != elements.end(); ++it)
996 {
997 size_t n = 0;
998 MPI_Unpack((void *) buf, bs, &p, &(n), 1, MPI_SIZE_T, MPI_COMM_WORLD);
999 for (size_t i = 0; i < n; ++i)
1000 {
1001 size_t index = 0;
1002 MPI_Unpack((void *) buf, bs, &p, &(index), 1, MPI_SIZE_T, MPI_COMM_WORLD);
1003 SymFncStatistics::Container& d = it->statistics.data[index];
1004 size_t countEW = 0;
1005 MPI_Unpack((void *) buf, bs, &p, &(countEW ), 1, MPI_SIZE_T, MPI_COMM_WORLD);
1006 MPI_Unpack((void *) buf, bs, &p, &(d.type ), 1, MPI_SIZE_T, MPI_COMM_WORLD);
1007 MPI_Unpack((void *) buf, bs, &p, &(d.Gmin ), 1, MPI_DOUBLE, MPI_COMM_WORLD);
1008 MPI_Unpack((void *) buf, bs, &p, &(d.Gmax ), 1, MPI_DOUBLE, MPI_COMM_WORLD);
1009 // d.element
1010 size_t ts = 0;
1011 MPI_Unpack((void *) buf, bs, &p, &ts , 1, MPI_SIZE_T, MPI_COMM_WORLD);
1012 char* element = new char[ts];
1013 MPI_Unpack((void *) buf, bs, &p, element , ts, MPI_CHAR , MPI_COMM_WORLD);
1014 d.element = element;
1015 delete[] element;
1016 // indexStructureEW.
1017 d.indexStructureEW.resize(d.countEW + countEW);
1018 MPI_Unpack((void *) buf, bs, &p, &(d.indexStructureEW[d.countEW]), countEW, MPI_SIZE_T, MPI_COMM_WORLD);
1019 // indexAtomEW.
1020 d.indexAtomEW.resize(d.countEW + countEW);
1021 MPI_Unpack((void *) buf, bs, &p, &(d.indexAtomEW[d.countEW] ), countEW, MPI_SIZE_T, MPI_COMM_WORLD);
1022 // valueEW.
1023 d.valueEW.resize(d.countEW + countEW);
1024 MPI_Unpack((void *) buf, bs, &p, &(d.valueEW[d.countEW] ), countEW, MPI_DOUBLE, MPI_COMM_WORLD);
1025
1026 d.countEW += countEW;
1027 }
1028 }
1029#endif
1030 return;
1031}
double d

References d, nnp::Mode::elements, MPI_SIZE_T, and p.

◆ writeExtrapolationWarnings()

void InterfaceLammps::writeExtrapolationWarnings ( )

Write extrapolation warnings to log.

Definition at line 1033 of file InterfaceLammps.cpp.

1034{
1035 for (vector<Element>::const_iterator it = elements.begin();
1036 it != elements.end(); ++it)
1037 {
1038 vector<string> vs = it->statistics.getExtrapolationWarningLines();
1039 for (vector<string>::const_iterator it2 = vs.begin();
1040 it2 != vs.end(); ++it2)
1041 {
1042 log << (*it2);
1043 }
1044 }
1045
1046 return;
1047}

References nnp::Mode::elements, and nnp::Mode::log.

Referenced by initialize().

Here is the caller graph for this function:

◆ clearExtrapolationWarnings()

void InterfaceLammps::clearExtrapolationWarnings ( )

Clear extrapolation warnings storage.

Definition at line 1049 of file InterfaceLammps.cpp.

1050{
1051 for (vector<Element>::iterator it = elements.begin();
1052 it != elements.end(); ++it)
1053 {
1054 it->statistics.clear();
1055 }
1056
1057 return;
1058}

References nnp::Mode::elements.

◆ addCharge()

void InterfaceLammps::addCharge ( int index,
double Q )

Read atomic charges from LAMMPS into n2p2.

Definition at line 648 of file InterfaceLammps.cpp.

649{
650 Atom& a = structure.atoms.at(index);
651 a.charge = Q;
652 //log << strpr("Atom %5zu (%2s) q: %24.16E\n",
653 // a.tag, elementMap[a.element].c_str(), a.charge);
654}
double charge
Atomic charge determined by neural network.
Definition Atom.h:121

References nnp::Atom::charge, and structure.

◆ getQEqParams()

void InterfaceLammps::getQEqParams ( double *const & atomChi,
double *const & atomJ,
double *const & sigmaSqrtPi,
double *const *const & gammaSqrt2,
double & qRef ) const

Write QEq arrays from n2p2 to LAMMPS.

Parameters
[in]atomChiElectronegativities.
[in]atomJAtomic hardness.
[in]atomSigmaGaussian width.
[in]qRefReference charge of the structure.

Definition at line 572 of file InterfaceLammps.cpp.

574{
575 for (size_t i = 0; i < structure.atoms.size(); ++i) {
576 Atom const& a = structure.atoms.at(i);
577 size_t const ia = a.index;
578 atomChi[ia] = a.chi;
579 }
580
581 for (size_t i = 0; i < numElements; ++i)
582 {
583 double const iSigma = elements.at(i).getQsigma();
584 atomJ[i] = elements.at(i).getHardness();
585 sigmaSqrtPi[i] = sqrt(M_PI) * iSigma;
586 for (size_t j = 0; j < numElements; j++)
587 {
588 double const jSigma = elements.at(j).getQsigma();
589 gammaSqrt2[i][j] = sqrt(2.0 * (iSigma * iSigma + jSigma * jSigma));
590 }
591 }
592 qRef = structure.chargeRef;
593}
double chi
Atomic electronegativity determined by neural network.
Definition Atom.h:119

References nnp::Atom::chi, nnp::Mode::elements, nnp::Atom::index, nnp::Mode::numElements, and structure.

◆ getdEdQ()

void InterfaceLammps::getdEdQ ( double *const & dEtotdQ) const

Write the derivative of total energy with respect to atomic charges from n2p2 into LAMMPS.

Parameters
[in]dEtotdQDerivative of the total energy w.r.t. atomic charge.

Definition at line 595 of file InterfaceLammps.cpp.

596{
597 Atom const* ai = NULL;
598 for (size_t i = 0; i < structure.atoms.size(); ++i)
599 {
600 ai = &(structure.atoms.at(i));
601 size_t const ia = ai->index;
602 dEtotdQ[ia] += ai->dEdG.back();
603 }
604}

References nnp::Atom::dEdG, nnp::Atom::index, and structure.

◆ getScreeningInfo()

void InterfaceLammps::getScreeningInfo ( double *const & rScreen) const

Read screening function information from n2p2 into LAMMPS.

Parameters
[in]rScreenArray that contains screening radii.

Definition at line 656 of file InterfaceLammps.cpp.

657{
658 screenInfo[0] = (double) screeningFunction.getCoreFunctionType(); //TODO: this does not work atm
659 screenInfo[1] = screeningFunction.getInner();
660 screenInfo[2] = screeningFunction.getOuter();
661 screenInfo[3] = 1.0 / (screenInfo[2] - screenInfo[1]); // scale
662}

References nnp::Mode::screeningFunction.

◆ getdChidxyz()

void InterfaceLammps::getdChidxyz ( int tag,
double *const & dChidx,
double *const & dChidy,
double *const & dChidz ) const

Transfer spatial derivatives of atomic electronegativities.

Parameters
[in]tagAtom of interest
dChidx
dChidy
dChidz

Definition at line 606 of file InterfaceLammps.cpp.

608{
609 Atom const &ai = structure.atoms.at(ind);
610
611 for (size_t j = 0; j < structure.numAtoms; ++j) {
612 Atom const &aj = structure.atoms.at(j);
613#ifndef NNP_FULL_SFD_MEMORY
614 vector <vector<size_t>> const &tableFull
615 = elements.at(aj.element).getSymmetryFunctionTable();
616#endif
617 Vec3D dChi;
618 // need to add this case because the loop over the neighbors
619 // does not include the contribution dChi_i/dr_i.
620 if (ai.index == j) {
621 for (size_t k = 0; k < aj.numSymmetryFunctions; ++k) {
622 dChi += aj.dChidG.at(k) * aj.dGdr.at(k);
623 }
624 }
625
626 for (auto const &n : aj.neighbors) {
627 if (n.d > maxCutoffRadius) break;
628 if (n.index == ai.index) {
629#ifndef NNP_FULL_SFD_MEMORY
630 vector <size_t> const &table = tableFull.at(n.element);
631 for (size_t k = 0; k < n.dGdr.size(); ++k) {
632 dChi += aj.dChidG.at(table.at(k)) * n.dGdr.at(k);
633 }
634#else
635 for (size_t k = 0; k < aj.numSymmetryFunctions; ++k)
636 {
637 dChi += aj.dChidG.at(k) * n.dGdr.at(k);
638 }
639#endif
640 }
641 }
642 dChidx[j] = dChi[0];
643 dChidy[j] = dChi[1];
644 dChidz[j] = dChi[2];
645 }
646}

References nnp::Atom::dChidG, nnp::Atom::dGdr, nnp::Atom::element, nnp::Mode::elements, nnp::Atom::index, nnp::Mode::maxCutoffRadius, nnp::Atom::neighbors, nnp::Atom::numSymmetryFunctions, and structure.

◆ setElecDone()

void InterfaceLammps::setElecDone ( )

Set isElecDone true after running the first NN in 4G-HDNNPs.

Definition at line 669 of file InterfaceLammps.cpp.

670{
671 if (isElecDone) isElecDone = false;
672}

References isElecDone.

◆ writeToFile()

void InterfaceLammps::writeToFile ( std::string const fileName,
bool const append )

Write current structure to file in units used in training data.

Parameters
fileNameFile name of the output structure file.
appendtrue if structure should be appended to existing file.

Definition at line 1060 of file InterfaceLammps.cpp.

1062{
1064 structure.writeToFile(fileName, false, append);
1065 structure.toNormalizedUnits(meanEnergy, convEnergy, convLength, convCharge);
1066}
double convCharge
Definition Mode.h:675

References nnp::Mode::convCharge, nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Mode::meanEnergy, and structure.

◆ add3DVecToArray()

void InterfaceLammps::add3DVecToArray ( double *const & arr,
Vec3D const & v ) const

Add a Vec3D vector to a 3D array in place.

Parameters
[in,out]arrArray which is edited in place.
[in]vVector which is added to arr.

Definition at line 1068 of file InterfaceLammps.cpp.

1069{
1070 arr[0] += v[0];
1071 arr[1] += v[1];
1072 arr[2] += v[2];
1073}

Referenced by getForcesDevelop().

Here is the caller graph for this function:

Member Data Documentation

◆ myRank

int nnp::InterfaceLammps::myRank
protected

Process rank.

Definition at line 289 of file InterfaceLammps.h.

Referenced by initialize(), InterfaceLammps(), and setLocalAtoms().

◆ initialized

bool nnp::InterfaceLammps::initialized
protected

Initialization state.

Definition at line 291 of file InterfaceLammps.h.

Referenced by initialize(), InterfaceLammps(), and isInitialized().

◆ hasGlobalStructure

bool nnp::InterfaceLammps::hasGlobalStructure
protected

Whether n2p2 knows about the global structure or only a local part.

Definition at line 293 of file InterfaceLammps.h.

Referenced by getGlobalStructureStatus(), InterfaceLammps(), and setGlobalStructureStatus().

◆ showew

bool nnp::InterfaceLammps::showew
protected

Corresponds to LAMMPS showew keyword.

Definition at line 295 of file InterfaceLammps.h.

Referenced by initialize(), and InterfaceLammps().

◆ resetew

bool nnp::InterfaceLammps::resetew
protected

Corresponds to LAMMPS resetew keyword.

Definition at line 297 of file InterfaceLammps.h.

Referenced by initialize(), and InterfaceLammps().

◆ showewsum

int nnp::InterfaceLammps::showewsum
protected

Corresponds to LAMMPS showewsum keyword.

Definition at line 299 of file InterfaceLammps.h.

Referenced by initialize(), and InterfaceLammps().

◆ maxew

int nnp::InterfaceLammps::maxew
protected

Corresponds to LAMMPS maxew keyword.

Definition at line 301 of file InterfaceLammps.h.

Referenced by initialize(), and InterfaceLammps().

◆ cflength

double nnp::InterfaceLammps::cflength
protected

◆ cfenergy

double nnp::InterfaceLammps::cfenergy
protected

Corresponds to LAMMPS cfenergy keyword.

Definition at line 305 of file InterfaceLammps.h.

Referenced by getAtomicEnergy(), getEnergy(), getForces(), getForcesChi(), getForcesDevelop(), initialize(), and InterfaceLammps().

◆ emap

std::string nnp::InterfaceLammps::emap
protected

Corresponds to LAMMPS map keyword.

Definition at line 307 of file InterfaceLammps.h.

Referenced by initialize().

◆ indexMap

std::vector<size_t> nnp::InterfaceLammps::indexMap
protected

Map from LAMMPS index to n2p2 atom index.

Definition at line 309 of file InterfaceLammps.h.

Referenced by addNeighbor(), and setLocalAtoms().

◆ ignoreType

std::map<int, bool> nnp::InterfaceLammps::ignoreType
protected

True if atoms of this LAMMPS type will be ignored.

Definition at line 311 of file InterfaceLammps.h.

Referenced by addNeighbor(), initialize(), and setLocalAtoms().

◆ mapTypeToElement

std::map<int, std::size_t> nnp::InterfaceLammps::mapTypeToElement
protected

Map from LAMMPS type to n2p2 element index.

Definition at line 313 of file InterfaceLammps.h.

Referenced by addNeighbor(), initialize(), and setLocalAtoms().

◆ mapElementToType

std::map<std::size_t, int> nnp::InterfaceLammps::mapElementToType
protected

Map from n2p2 element index to LAMMPS type.

Definition at line 315 of file InterfaceLammps.h.

Referenced by initialize().

◆ structure

◆ isElecDone

bool nnp::InterfaceLammps::isElecDone
protected

True if first NN is calculated.

Definition at line 319 of file InterfaceLammps.h.

Referenced by InterfaceLammps(), process(), and setElecDone().


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