n2p2 - A neural network potential package
Loading...
Searching...
No Matches
Structure.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 STRUCTURE_H
18#define STRUCTURE_H
19
20#include "Atom.h"
21#include "Element.h"
22#include "ElementMap.h"
23#include "ErfcBuf.h"
24#include "EwaldSetup.h"
25#include "ScreeningFunction.h"
26#include "Vec3D.h"
27#include <Eigen/Core> // MatrixXd, VectorXd
28#include <cstddef> // std::size_t
29#include <fstream> // std::ofstream
30#include <map> // std::map
31#include <string> // std::string
32#include <vector> // std::vector
33
34namespace nnp
35{
36
39{
57
73 std::size_t index;
75 std::size_t numAtoms;
77 std::size_t numElements;
79 std::size_t numElementsPresent;
81 int pbc[3];
83 double energy;
85 double energyRef;
86 // TODO: MPI_Pack
89 // TODO: MPI_Pack
91 double energyElec;
96 double charge;
98 double chargeRef;
100 double volume;
104 // TODO: MPI_Pack
106 double lambda;
110 std::string comment;
116 Eigen::MatrixXd A;
120 std::vector<std::size_t> numAtomsPerElement;
122 std::vector<Atom> atoms;
123
126 Structure();
143 void addAtom(Atom const& atom,
144 std::string const& element);
151 void readFromFile(std::string const fileName
152 = "input.data");
160 void readFromFile(std::ifstream& file);
167 void readFromLines(std::vector<
168 std::string> const& lines);
178 EwaldSetup& ewaldSetup,
179 double rcutScreen,
180 double maxCutoffRadius);
189 double cutoffRadius,
190 bool sortByDistance = false);
200 double cutoffRadius,
201 std::vector<
202 std::vector<double>>& cutoffs);
205 void sortNeighborList();
212 void setupNeighborCutoffMap(std::vector<
213 std::vector<double>> cutoffs);
221 void calculatePbcCopies(double cutoffRadius, int (&pbc)[3]);
262 void calculateInverseBox();
269 bool canMinimumImageConventionBeApplied(double cutoffRadius);
278 void calculateVolume();
294 EwaldSetup& ewaldSetup,
295 Eigen::VectorXd hardness,
296 Eigen::MatrixXd gammaSqrt2,
297 Eigen::VectorXd sigmaSqrtPi,
298 ScreeningFunction const& fs,
299 double const fourPiEps,
300 ErfcBuf& erfcBuf);
315 Eigen::MatrixXd gammaSqrt2,
316 Eigen::VectorXd sigmaSqrtPi,
317 ScreeningFunction const& fs,
318 double const fourPiEps);
329 void calculateDAdrQ(
330 EwaldSetup& ewaldSetup,
331 Eigen::MatrixXd gammaSqrt2,
332 double const fourPiEps,
333 ErfcBuf& erfcBuf);
338 void calculateDQdChi(std::vector<Eigen::VectorXd> &dQdChi);
343 void calculateDQdJ(std::vector<Eigen::VectorXd> &dQdJ);
354 void calculateDQdr(
355 std::vector<size_t> const& atomsIndices,
356 std::vector<size_t> const& compIndices,
357 double const maxCutoffRadius,
358 std::vector<Element> const& elements);
372 Eigen::VectorXd hardness,
373 Eigen::MatrixXd gammaSqrt2,
374 Eigen::VectorXd sigmaSqrtPi,
375 ScreeningFunction const& fs,
376 double const fourPiEps);
380 Eigen::VectorXd const calculateForceLambdaTotal() const;
384 Eigen::VectorXd const calculateForceLambdaElec() const;
387 void remap();
392 void remap(Atom& atom);
401 void toNormalizedUnits(double meanEnergy,
402 double convEnergy,
403 double convLength,
404 double convCharge);
413 void toPhysicalUnits(double meanEnergy,
414 double convEnergy,
415 double convLength,
416 double convCharge);
421 std::size_t getMaxNumNeighbors() const;
427 void freeAtoms(bool all,
428 double const maxCutoffRadius = 0.0);
431 void reset();
434 void clearNeighborList();
438 void clearElectrostatics(bool clearDQdr = false);
465 void updateError(
466 std::string const& property,
467 std::map<std::string, double>& error,
468 std::size_t& count) const;
473 std::string getEnergyLine() const;
478 std::vector<std::string> getForcesLines() const;
483 std::vector<std::string> getChargesLines() const;
491 void writeToFile(
492 std::string const fileName ="output.data",
493 bool const ref = true,
494 bool const append = false) const;
501 void writeToFile(
502 std::ofstream* const& file,
503 bool const ref = true) const;
508 void writeToFileXyz(std::ofstream* const& file) const;
517 std::ofstream* const& file) const;
524 std::ofstream* const& file,
525 std::string const elements) const;
530 std::vector<std::string> info() const;
531};
532
533}
534
535#endif
Contains element map.
Definition ElementMap.h:30
Setup data for Ewald summation.
Definition EwaldSetup.h:37
A screening functions for use with electrostatics.
Definition Atom.h:29
Storage for a single atom.
Definition Atom.h:33
Helper class to store previously calculated values of erfc() that are needed during the charge equili...
Definition ErfcBuf.h:17
void freeAtoms(bool all, double const maxCutoffRadius=0.0)
Free symmetry function memory for all atoms, see free() in Atom class.
void calculateVolume()
Calculate volume from box vectors.
SampleType
Enumerates different sample types (e.g.
Definition Structure.h:43
@ ST_VALIDATION
Structure is part of validation set (currently unused).
Definition Structure.h:52
@ ST_TRAINING
Structure is part of the training set.
Definition Structure.h:49
@ ST_UNKNOWN
Sample type not assigned yet.
Definition Structure.h:46
@ ST_TEST
Structure is part of the test set.
Definition Structure.h:55
Vec3D invbox[3]
Inverse simulation box vectors.
Definition Structure.h:114
void toPhysicalUnits(double meanEnergy, double convEnergy, double convLength, double convCharge)
Switch to physical units, shift energy and change energy, length and charge unit.
Vec3D box[3]
Simulation box vectors.
Definition Structure.h:112
bool isTriclinic
If the simulation box is triclinic.
Definition Structure.h:63
std::size_t getMaxNumNeighbors() const
Find maximum number of neighbors.
void setupNeighborCutoffMap(std::vector< std::vector< double > > cutoffs)
Set up a neighbor cut-off map which gives the index (value) of the last needed neighbor corresponding...
void sortNeighborList()
Sort all neighbor lists of this structure with respect to the distance.
Eigen::VectorXd const calculateForceLambdaElec() const
Calculate lambda_elec vector which is needed for the electrostatic force calculation in 4G NN.
std::vector< std::size_t > numAtomsPerElement
Number of atoms of each element in this structure.
Definition Structure.h:120
double lambda
Lagrange multiplier used for charge equilibration.
Definition Structure.h:106
void calculateMaxCutoffRadiusOverall(EwaldSetup &ewaldSetup, double rcutScreen, double maxCutoffRadius)
Calculate maximal cut-off if cut-off of screening and real part Ewald summation are also considered.
bool canMinimumImageConventionBeApplied(double cutoffRadius)
Check if cut-off radius is small enough to apply minimum image convention.
void writeToFile(std::string const fileName="output.data", bool const ref=true, bool const append=false) const
Write configuration to file.
void calculateDQdChi(std::vector< Eigen::VectorXd > &dQdChi)
Calculates derivative of the charges with respect to electronegativities.
std::string comment
Structure comment.
Definition Structure.h:110
bool isPeriodic
If periodic boundary conditions apply.
Definition Structure.h:61
double maxCutoffRadiusOverall
Maximum cut-off radius with respect to symmetry functions, screening function and Ewald summation.
Definition Structure.h:103
double charge
Charge determined by neural network potential.
Definition Structure.h:96
double calculateElectrostaticEnergy(EwaldSetup &ewaldSetup, Eigen::VectorXd hardness, Eigen::MatrixXd gammaSqrt2, Eigen::VectorXd sigmaSqrtPi, ScreeningFunction const &fs, double const fourPiEps, ErfcBuf &erfcBuf)
Compute electrostatic energy with global charge equilibration.
bool NeighborListIsSorted
If the neighbor list has been sorted by distance.
Definition Structure.h:67
ElementMap elementMap
Copy of element map provided as constructor argument.
Definition Structure.h:59
void setElementMap(ElementMap const &elementMap)
Set element map of structure.
Definition Structure.cpp:71
std::size_t index
Index number of this structure.
Definition Structure.h:73
void remap()
Translate all atoms back into box if outside.
std::vector< std::string > getChargesLines() const
Get reference and NN charges for all atoms.
double chargeRef
Reference charge.
Definition Structure.h:98
SampleType sampleType
Sample type (training or test set).
Definition Structure.h:108
double energyShort
Short-range part of the potential energy predicted by NNP.
Definition Structure.h:88
void calculateDQdJ(std::vector< Eigen::VectorXd > &dQdJ)
Calculates derivative of the charges with respect to atomic hardness.
bool hasAMatrix
If A matrix of this structure is currently stored.
Definition Structure.h:118
void addAtom(Atom const &atom, std::string const &element)
Add a single atom to structure.
Definition Structure.cpp:80
void calculateInverseBox()
Calculate inverse box.
void writeToFilePoscar(std::ofstream *const &file) const
Write configuration to POSCAR file.
void writeToFileXyz(std::ofstream *const &file) const
Write configuration to xyz file.
bool hasSymmetryFunctionDerivatives
If symmetry function derivatives are saved for each atom.
Definition Structure.h:71
void calculateDQdr(std::vector< size_t > const &atomsIndices, std::vector< size_t > const &compIndices, double const maxCutoffRadius, std::vector< Element > const &elements)
Calculates derivative of the charges with respect to the atom's position.
Eigen::MatrixXd A
Global charge equilibration matrix A'.
Definition Structure.h:116
void calculatePbcCopies(double cutoffRadius, int(&pbc)[3])
Calculate required PBC copies.
Eigen::VectorXd const calculateForceLambdaTotal() const
Calculate lambda_total vector which is needed for the total force calculation in 4G NN.
std::string getEnergyLine() const
Get reference and NN energy.
void clearElectrostatics(bool clearDQdr=false)
Clear A-matrix, dAdrQ and optionally dQdr.
double energyElec
Electrostatics part of the potential energy predicted by NNP.
Definition Structure.h:91
double energy
Potential energy determined by neural network.
Definition Structure.h:83
void readFromFile(std::string const fileName="input.data")
Read configuration from file.
Definition Structure.cpp:97
double energyRef
Reference potential energy.
Definition Structure.h:85
int pbc[3]
Number of PBC images necessary in each direction for max cut-off.
Definition Structure.h:81
void toNormalizedUnits(double meanEnergy, double convEnergy, double convLength, double convCharge)
Normalize structure, shift energy and change energy, length and charge unit.
void reset()
Reset everything but elementMap.
std::vector< std::string > getForcesLines() const
Get reference and NN forces for all atoms.
Structure()
Constructor, initializes to zero.
Definition Structure.cpp:35
double volume
Simulation box volume.
Definition Structure.h:100
void calculateElectrostaticEnergyDerivatives(Eigen::VectorXd hardness, Eigen::MatrixXd gammaSqrt2, Eigen::VectorXd sigmaSqrtPi, ScreeningFunction const &fs, double const fourPiEps)
Calculates partial derivatives of electrostatic Energy with respect to atom's coordinates and charges...
Vec3D applyMinimumImageConvention(Vec3D const &dr)
Calculate distance between two atoms in the minimum image convention.
void calculateNeighborList(double cutoffRadius, bool sortByDistance=false)
Calculate neighbor list for all atoms.
std::size_t numAtoms
Total number of atoms present in this structure.
Definition Structure.h:75
std::size_t numElements
Global number of elements (all structures).
Definition Structure.h:77
void updateError(std::string const &property, std::map< std::string, double > &error, std::size_t &count) const
Update property error metrics with this structure.
void calculateDAdrQ(EwaldSetup &ewaldSetup, Eigen::MatrixXd gammaSqrt2, double const fourPiEps, ErfcBuf &erfcBuf)
Calculates derivative of A-matrix with respect to the atoms positions and contract it with Q.
std::vector< Atom > atoms
Vector of all atoms in this structure.
Definition Structure.h:122
std::size_t numElementsPresent
Number of elements present in this structure.
Definition Structure.h:79
void readFromLines(std::vector< std::string > const &lines)
Read configuration from lines.
double calculateScreeningEnergy(Eigen::MatrixXd gammaSqrt2, Eigen::VectorXd sigmaSqrtPi, ScreeningFunction const &fs, double const fourPiEps)
Calculate screening energy which needs to be added (!) to the electrostatic energy in order to remove...
bool hasNeighborList
If the neighbor list has been calculated.
Definition Structure.h:65
void clearNeighborList()
Clear neighbor list of all atoms.
std::vector< std::string > info() const
Get structure information as a vector of strings.
bool hasCharges
If all charges of this structure have been calculated (and stay the same, e.g.
Definition Structure.h:94
bool hasSymmetryFunctions
If symmetry function values are saved for each atom.
Definition Structure.h:69
Vector in 3 dimensional real space.
Definition Vec3D.h:30