n2p2 - A neural network potential package
nnp-norm2.cpp File Reference
#include "Training.h"
#include "mpi-extra.h"
#include "utility.h"
#include <mpi.h>
#include <cmath>
#include <cstdlib>
#include <iostream>
#include <fstream>
#include <vector>
Include dependency graph for nnp-norm2.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 30 of file nnp-norm2.cpp.

31{
32 int numProcs = 0;
33 int myRank = 0;
34 size_t stage = 0;
35 ofstream myLog;
36
37 if (argc > 2)
38 {
39 cout << "USAGE: " << argv[0] << " <stage>\n"
40 << " <stage> ... Training stage (only required if training"
41 " is a multi-stage process.\n"
42 << " Execute in directory with these NNP files present:\n"
43 << " - input.data (structure file)\n"
44 << " - input.nn (NNP settings)\n"
45 << " - scaling.data (symmetry function scaling data)\n"
46 << " - \"weights(e).%%03d.data\" (weights files)\n";
47 return 1;
48 }
49
50 string suffix = "";
51 if (argc > 1)
52 {
53 stage = (size_t)atoi(argv[1]);
54 suffix = strpr(".stage-%zu", stage);
55 }
56
57 MPI_Init(&argc, &argv);
58 MPI_Comm_size(MPI_COMM_WORLD, &numProcs);
59 MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
60
61 Training training;
62 if (myRank != 0) training.log.writeToStdout = false;
63 myLog.open((strpr("nnp-norm2.log.%04d", myRank) + suffix).c_str());
64 training.log.registerStreamPointer(&myLog);
65 training.setupMPI();
66 training.initialize();
67 training.loadSettingsFile();
68 training.setStage(stage);
69 training.setupGeneric();
71 training.setupSymmetryFunctionStatistics(false, false, false, false);
72
73 auto nnpType = training.getNnpType();
74 if (nnpType == Training::NNPType::SHORT_CHARGE_NN && stage == 1)
75 {
76 throw runtime_error("ERROR: Normalization of charges not yet "
77 "implemented\n.");
78 }
79
80 if (training.settingsKeywordExists("mean_energy") ||
81 training.settingsKeywordExists("conv_energy") ||
82 training.settingsKeywordExists("conv_length"))
83 {
84 throw runtime_error("ERROR: Normalization keywords found in settings, "
85 "please remove them first.\n");
86 }
87
88 if (!training.settingsKeywordExists("use_short_forces"))
89 {
90 throw runtime_error("ERROR: Normalization is only possible if forces "
91 "are used (keyword \"use_short_forces\").\n");
92 }
93
94 // Need RNG for initial random weights.
96
97 // Distribute structures to MPI processes.
98 training.distributeStructures(false);
99
100 // Initialize weights and biases for neural networks.
101 training.initializeWeights();
102
103 training.writeWeights("short", "weights.%03zu.norm");
104 if (nnpType == Training::NNPType::SHORT_CHARGE_NN && stage == 2)
105 {
106 training.writeWeights("charge", "weightse.%03zu.norm");
107 }
108
109 training.log << "\n";
110 training.log << "*** DATA SET NORMALIZATION **************"
111 "**************************************\n";
112 training.log << "\n";
113
114 ofstream fileEvsV;
115 fileEvsV.open(strpr("evsv.dat.%04d", myRank).c_str());
116 if (myRank == 0)
117 {
118 // File header.
119 vector<string> title;
120 vector<string> colName;
121 vector<string> colInfo;
122 vector<size_t> colSize;
123 title.push_back("Energy vs. volume comparison.");
124 colSize.push_back(16);
125 colName.push_back("V_atom");
126 colInfo.push_back("Volume per atom.");
127 colSize.push_back(16);
128 colName.push_back("Eref_atom");
129 colInfo.push_back("Reference energy per atom.");
130 colSize.push_back(10);
131 colName.push_back("N");
132 colInfo.push_back("Number of atoms.");
133 colSize.push_back(16);
134 colName.push_back("V");
135 colInfo.push_back("Volume of structure.");
136 colSize.push_back(16);
137 colName.push_back("Eref");
138 colInfo.push_back("Reference energy of structure.");
139 colSize.push_back(16);
140 colName.push_back("Eref_offset");
141 colInfo.push_back("Reference energy of structure (including offset).");
142 appendLinesToFile(fileEvsV,
143 createFileHeader(title, colSize, colName, colInfo));
144 }
145
146 size_t numAtomsTotal = 0;
147 size_t numStructures = 0;
148 double meanEnergyPerAtomRef = 0.0;
149 double meanEnergyPerAtomNnp = 0.0;
150 double sigmaEnergyPerAtomRef = 0.0;
151 double sigmaEnergyPerAtomNnp = 0.0;
152 double meanForceRef = 0.0;
153 double meanForceNnp = 0.0;
154 double sigmaForceRef = 0.0;
155 double sigmaForceNnp = 0.0;
156 training.log << "Computing initial prediction for all structures...\n";
157 for (auto& s : training.structures)
158 {
159 // File output for evsv.dat.
160 fileEvsV << strpr("%16.8E %16.8E %10zu %16.8E %16.8E %16.8E\n",
161 s.volume / s.numAtoms,
162 s.energyRef / s.numAtoms,
163 s.numAtoms,
164 s.volume,
165 s.energyRef,
166 training.getEnergyWithOffset(s, true));
167 s.calculateNeighborList(training.getMaxCutoffRadius());
168#ifdef N2P2_NO_SF_GROUPS
169 training.calculateSymmetryFunctions(s, true);
170#else
171 training.calculateSymmetryFunctionGroups(s, true);
172#endif
173 training.calculateAtomicNeuralNetworks(s, true);
174 training.calculateEnergy(s);
175 training.calculateForces(s);
176 s.clearNeighborList();
177
178 numStructures++;
179 numAtomsTotal += s.numAtoms;
180 meanEnergyPerAtomRef += s.energyRef / s.numAtoms;
181 meanEnergyPerAtomNnp += s.energy / s.numAtoms;
182 for (auto& a : s.atoms)
183 {
184 meanForceRef += a.fRef[0] + a.fRef[1] + a.fRef[2];
185 meanForceNnp += a.f [0] + a.f [1] + a.f [2];
186 }
187 }
188
189 fileEvsV.flush();
190 fileEvsV.close();
191 MPI_Barrier(MPI_COMM_WORLD);
192 training.log << "Writing energy/atom vs. volume/atom data "
193 << "to \"evsv.dat\".\n";
194 if (myRank == 0) training.combineFiles("evsv.dat");
195 MPI_Allreduce(MPI_IN_PLACE, &numStructures , 1, MPI_SIZE_T, MPI_SUM, MPI_COMM_WORLD);
196 MPI_Allreduce(MPI_IN_PLACE, &numAtomsTotal , 1, MPI_SIZE_T, MPI_SUM, MPI_COMM_WORLD);
197 MPI_Allreduce(MPI_IN_PLACE, &meanEnergyPerAtomRef, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
198 MPI_Allreduce(MPI_IN_PLACE, &meanEnergyPerAtomNnp, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
199 MPI_Allreduce(MPI_IN_PLACE, &meanForceRef , 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
200 MPI_Allreduce(MPI_IN_PLACE, &meanForceNnp , 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
201 meanEnergyPerAtomRef /= numStructures;
202 meanEnergyPerAtomNnp /= numStructures;
203 meanForceRef /= 3 * numAtomsTotal;
204 meanForceNnp /= 3 * numAtomsTotal;
205 for (auto const& s : training.structures)
206 {
207 double ediffRef = s.energyRef / s.numAtoms - meanEnergyPerAtomRef;
208 double ediffNnp = s.energy / s.numAtoms - meanEnergyPerAtomNnp;
209 sigmaEnergyPerAtomRef += ediffRef * ediffRef;
210 sigmaEnergyPerAtomNnp += ediffNnp * ediffNnp;
211 for (auto const& a : s.atoms)
212 {
213 double fdiffRef = a.fRef[0] - meanForceRef;
214 double fdiffNnp = a.f [0] - meanForceNnp;
215 sigmaForceRef += fdiffRef * fdiffRef;
216 sigmaForceNnp += fdiffNnp * fdiffNnp;
217 fdiffRef = a.fRef[1] - meanForceRef;
218 fdiffNnp = a.f [1] - meanForceNnp;
219 sigmaForceRef += fdiffRef * fdiffRef;
220 sigmaForceNnp += fdiffNnp * fdiffNnp;
221 fdiffRef = a.fRef[2] - meanForceRef;
222 fdiffNnp = a.f [2] - meanForceNnp;
223 sigmaForceRef += fdiffRef * fdiffRef;
224 sigmaForceNnp += fdiffNnp * fdiffNnp;
225 }
226 }
227 MPI_Allreduce(MPI_IN_PLACE, &sigmaEnergyPerAtomRef, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
228 MPI_Allreduce(MPI_IN_PLACE, &sigmaEnergyPerAtomNnp, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
229 MPI_Allreduce(MPI_IN_PLACE, &sigmaForceRef , 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
230 MPI_Allreduce(MPI_IN_PLACE, &sigmaForceNnp , 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
231 sigmaEnergyPerAtomRef = sqrt(sigmaEnergyPerAtomRef / (numStructures - 1));
232 sigmaEnergyPerAtomNnp = sqrt(sigmaEnergyPerAtomNnp / (numStructures - 1));
233 sigmaForceRef = sqrt(sigmaForceRef / (3 * numAtomsTotal - 1));
234 sigmaForceNnp = sqrt(sigmaForceNnp / (3 * numAtomsTotal - 1));
235 training.log << "\n";
236 training.log << strpr("Total number of structures : %zu\n", numStructures);
237 training.log << strpr("Total number of atoms : %zu\n", numAtomsTotal);
238 training.log << "----------------------------------\n";
239 training.log << "Reference data statistics:\n";
240 training.log << "----------------------------------\n";
241 training.log << strpr("Mean/sigma energy per atom : %16.8E +/- %16.8E\n",
242 meanEnergyPerAtomRef,
243 sigmaEnergyPerAtomRef);
244 training.log << strpr("Mean/sigma force : %16.8E +/- %16.8E\n",
245 meanForceRef,
246 sigmaForceRef);
247 training.log << "----------------------------------\n";
248 training.log << "Initial NNP prediction statistics:\n";
249 training.log << "----------------------------------\n";
250 training.log << strpr("Mean/sigma energy per atom : %16.8E +/- %16.8E\n",
251 meanEnergyPerAtomNnp,
252 sigmaEnergyPerAtomNnp);
253 training.log << strpr("Mean/sigma force : %16.8E +/- %16.8E\n",
254 meanForceNnp,
255 sigmaForceNnp);
256 training.log << "----------------------------------\n";
257 double convEnergy = sigmaForceNnp / sigmaForceRef;
258 double convLength = sigmaForceNnp;
259 training.log << strpr("Conversion factor energy : %24.16E\n", convEnergy);
260 training.log << strpr("Conversion factor length : %24.16E\n", convLength);
261
262 //ofstream fileCfg;
263 //fileCfg.open(strpr("output.data.%04d", myRank).c_str());
264 //for (vector<Structure>::iterator it = training.structures.begin();
265 // it != training.structures.end(); ++it)
266 //{
267 // it->energyRef = (it->energyRef - meanEnergyPerAtom * it->numAtoms)
268 // * convEnergy;
269 // it->box[0] *= convLength;
270 // it->box[1] *= convLength;
271 // it->box[2] *= convLength;
272 // for (vector<Atom>::iterator it2 = it->atoms.begin();
273 // it2 != it->atoms.end(); ++it2)
274 // {
275 // it2->r *= convLength;
276 // it2->fRef *= convEnergy / convLength;
277 // }
278 // it->writeToFile(&fileCfg);
279 //}
280 //fileCfg.flush();
281 //fileCfg.close();
282 //MPI_Barrier(MPI_COMM_WORLD);
283 //training.log << "\n";
284 //training.log << "Writing converted data file to \"output.data\".\n";
285 //training.log << "WARNING: This data set is provided for debugging "
286 // "purposes only and is NOT intended for training.\n";
287 //if (myRank == 0) training.combineFiles("output.data");
288
289 if (myRank == 0)
290 {
291 training.log << "\n";
292 training.log << "Writing backup of original settings file to "
293 "\"input.nn.bak\".\n";
294 ofstream fileSettings;
295 fileSettings.open("input.nn.bak");
296 training.writeSettingsFile(&fileSettings);
297 fileSettings.close();
298
299 training.log << "\n";
300 training.log << "Writing extended settings file to \"input.nn\".\n";
301 training.log << "Use this settings file for normalized training.\n";
302 fileSettings.open("input.nn");
303 fileSettings << "#########################################"
304 "######################################\n";
305 fileSettings << "# DATA SET NORMALIZATION\n";
306 fileSettings << "#########################################"
307 "######################################\n";
308 fileSettings << strpr("mean_energy %24.16E # nnp-norm2\n",
309 meanEnergyPerAtomRef);
310 fileSettings << strpr("conv_energy %24.16E # nnp-norm2\n", convEnergy);
311 fileSettings << strpr("conv_length %24.16E # nnp-norm2\n", convLength);
312 fileSettings << "#########################################"
313 "######################################\n";
314 fileSettings << "\n";
315 training.writeSettingsFile(&fileSettings);
316 fileSettings.close();
317 }
318
319 training.log << "*****************************************"
320 "**************************************\n";
321
322 myLog.close();
323
324 MPI_Finalize();
325
326 return 0;
327}
int distributeStructures(bool randomize, bool excludeRank0=false, std::string const &fileName="input.data")
Read data file and distribute structures among processors.
Definition: Dataset.cpp:722
void setupMPI()
Initialize MPI with MPI_COMM_WORLD.
Definition: Dataset.cpp:52
std::vector< Structure > structures
All structures in this dataset.
Definition: Dataset.h:195
void setupRandomNumberGenerator()
Initialize random number generator.
Definition: Dataset.cpp:110
void combineFiles(std::string filePrefix) const
Combine individual MPI proc files to one.
Definition: Dataset.cpp:1742
void registerStreamPointer(std::ofstream *const &streamPointer)
Register new C++ ofstream pointer.
Definition: Log.cpp:91
bool writeToStdout
Turn on/off output to stdout.
Definition: Log.h:85
NNPType getNnpType() const
Getter for Mode::nnpType.
Definition: Mode.h:543
void initialize()
Write welcome message with version information.
Definition: Mode.cpp:50
double getMaxCutoffRadius() const
Getter for Mode::maxCutoffRadius.
Definition: Mode.h:563
void calculateAtomicNeuralNetworks(Structure &structure, bool const derivatives)
Calculate a single atomic neural network for a given atom and nn type.
Definition: Mode.cpp:1310
void calculateEnergy(Structure &structure) const
Calculate potential energy for a given structure.
Definition: Mode.cpp:1366
void calculateSymmetryFunctionGroups(Structure &structure, bool const derivatives)
Calculate all symmetry function groups for all atoms in given structure.
Definition: Mode.cpp:1230
virtual void setupSymmetryFunctionScaling(std::string const &fileName="scaling.data")
Set up symmetry function scaling from file.
Definition: Mode.cpp:532
bool settingsKeywordExists(std::string const &keyword) const
Check if keyword was found in settings file.
Definition: Mode.cpp:1622
Log log
Global log file.
Definition: Mode.h:510
double getEnergyWithOffset(Structure const &structure, bool ref=true) const
Add atomic energy offsets and return energy.
Definition: Mode.cpp:1512
void setupGeneric(bool skipNormalize=false)
Combine multiple setup routines and provide a basic NNP setup.
Definition: Mode.cpp:197
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
Definition: Mode.cpp:1150
void calculateForces(Structure &structure) const
Calculate forces for all atoms in given structure.
Definition: Mode.cpp:1392
void setupSymmetryFunctionStatistics(bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
Set up symmetry function statistics collection.
Definition: Mode.cpp:923
void loadSettingsFile(std::string const &fileName="input.nn")
Open settings file and load all keywords into memory.
Definition: Mode.cpp:156
void writeSettingsFile(std::ofstream *const &file) const
Write complete settings file.
Definition: Mode.cpp:1650
Training methods.
Definition: Training.h:36
void initializeWeights()
Initialize weights for all elements.
Definition: Training.cpp:259
void writeWeights(std::string const &nnName, std::string const &fileNameFormat) const
Write weights to files (one file for each element).
Definition: Training.cpp:1559
void setStage(std::size_t stage)
Set training stage (if multiple stages are needed for NNP type).
Definition: Training.cpp:361
#define MPI_SIZE_T
Definition: mpi-extra.h:22
string strpr(const char *format,...)
String version of printf function.
Definition: utility.cpp:90
vector< string > createFileHeader(vector< string > const &title, vector< size_t > const &colSize, vector< string > const &colName, vector< string > const &colInfo, char const &commentChar)
Definition: utility.cpp:104
void appendLinesToFile(ofstream &file, vector< string > const lines)
Append multiple lines of strings to open file stream.
Definition: utility.cpp:218

References nnp::appendLinesToFile(), nnp::Mode::calculateAtomicNeuralNetworks(), nnp::Mode::calculateEnergy(), nnp::Mode::calculateForces(), nnp::Mode::calculateSymmetryFunctionGroups(), nnp::Mode::calculateSymmetryFunctions(), nnp::Dataset::combineFiles(), nnp::createFileHeader(), nnp::Dataset::distributeStructures(), nnp::Mode::getEnergyWithOffset(), nnp::Mode::getMaxCutoffRadius(), nnp::Mode::getNnpType(), nnp::Mode::initialize(), nnp::Training::initializeWeights(), nnp::Mode::loadSettingsFile(), nnp::Mode::log, MPI_SIZE_T, nnp::Log::registerStreamPointer(), nnp::Training::setStage(), nnp::Mode::settingsKeywordExists(), nnp::Mode::setupGeneric(), nnp::Dataset::setupMPI(), nnp::Dataset::setupRandomNumberGenerator(), nnp::Mode::setupSymmetryFunctionScaling(), nnp::Mode::setupSymmetryFunctionStatistics(), nnp::strpr(), nnp::Dataset::structures, nnp::Mode::writeSettingsFile(), nnp::Log::writeToStdout, and nnp::Training::writeWeights().

Here is the call graph for this function: