n2p2 - A neural network potential package
nnp-checkdw.cpp File Reference
#include "Training.h"
#include "mpi-extra.h"
#include "utility.h"
#include <mpi.h>
#include <algorithm>
#include <cmath>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <string>
#include <vector>
Include dependency graph for nnp-checkdw.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 33 of file nnp-checkdw.cpp.

34{
35 int warning = 0;
36 int numProcs = 0;
37 int myRank = 0;
38 size_t stage = 0;
39 double delta = 0.0;
40 ofstream myLog;
41 ofstream outFileSummary;
42 map<string, double> meanAbsError;
43 map<string, double> maxAbsError;
44 map<string, string> outFilesNames;
45 map<string, ofstream*> outFiles;
46
47 if (argc < 2 || argc > 3)
48 {
49 cout << "USAGE: " << argv[0] << " <stage> <<delta>>\n"
50 << " <stage> ..... Training stage (irrelevant for NNPs "
51 "with only one stage).\n"
52 << " <<delta>> ... (optional) Displacement for central "
53 "difference (default: 1.0e-4).\n"
54 << " Execute in directory with these NNP files present:\n"
55 << " - input.data (structure file)\n"
56 << " - input.nn (NNP settings)\n"
57 << " - scaling.data (symmetry function scaling data)\n"
58 << " - \"weights.%%03d.data\" (weights files)\n";
59 return 1;
60 }
61
62 MPI_Init(&argc, &argv);
63 MPI_Comm_size(MPI_COMM_WORLD, &numProcs);
64 MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
65
66 stage = (size_t)atoi(argv[1]);
67 if (argc == 3) delta = atof(argv[2]);
68 else delta = 1.0E-4;
69
70 Training training;
71 if (myRank != 0) training.log.writeToStdout = false;
72 myLog.open(strpr("nnp-checkdw.log.%04d", myRank).c_str());
73 training.log.registerStreamPointer(&myLog);
74 training.setupMPI();
75 training.initialize();
76 training.loadSettingsFile();
77 training.setStage(stage);
78 training.setupGeneric();
79 bool normalize = training.useNormalization();
81 training.setupSymmetryFunctionStatistics(false, false, true, false);
83 training.distributeStructures(false);
84 if (normalize) training.toNormalizedUnits();
85 auto pk = training.setupNumericDerivCheck();
86
87 training.log << "\n";
88 training.log << "*** ANALYTIC/NUMERIC WEIGHT DERIVATIVES C"
89 "HECK *********************************\n";
90 training.log << "\n";
91 training.log << strpr("Delta for symmetric difference quotient: %11.3E\n",
92 delta);
93
94 string fileName = "checkdw-summary.out";
95 training.log << strpr("Per-structure summary of analytic/numeric "
96 "weight derivative comparison\n"
97 "will be written to \"%s\"\n",
98 fileName.c_str());
99 fileName += strpr(".%04d", myRank);
100 outFileSummary.open(fileName.c_str());
101 if (myRank == 0)
102 {
103 // File header.
104 vector<string> title;
105 vector<string> colName;
106 vector<string> colInfo;
107 vector<size_t> colSize;
108 title.push_back(strpr("Per-structure summary of analytic vs. numeric "
109 "weight derivative comparison (delta = %11.3E).",
110 delta));
111 colSize.push_back(10);
112 colName.push_back("struct");
113 colInfo.push_back("Structure index.");
114 colSize.push_back(10);
115 colName.push_back("numAtoms");
116 colInfo.push_back("Number of atoms in structure.");
117 for (auto k : pk)
118 {
119 colSize.push_back(16);
120 colName.push_back("MAE_" + k);
121 colInfo.push_back("Mean over all absolute differences between "
122 "analytic and numeric weight derivatives for "
123 "property \"" + k + "\".");
124 colSize.push_back(16);
125 colName.push_back("maxAE_" + k);
126 colInfo.push_back("Maximum over all absolute differences between "
127 "analytic and numeric weight derivatives for "
128 "property \"" + k + "\".");
129 }
130 appendLinesToFile(outFileSummary,
131 createFileHeader(title, colSize, colName, colInfo));
132 }
133
134 for (auto k : pk)
135 {
136 string fileName = "checkdw-weights." + k + ".out";
137 training.log << "Individual analytic/numeric weight derivatives for "
138 "property \"" + k + "\"\nwill be written to \""
139 + fileName + "\"\n";
140 outFilesNames[k] = fileName;
141 fileName += strpr(".%04d", myRank);
142 outFiles[k] = new ofstream();
143 outFiles.at(k)->open(fileName.c_str());
144 if (myRank == 0)
145 {
146 // File header.
147 vector<string> title;
148 vector<string> colName;
149 vector<string> colInfo;
150 vector<size_t> colSize;
151 title.push_back(strpr("Comparison of analytic and numeric weight "
152 "derivatives for property \"%s\" (delta = "
153 "%11.3E).", k.c_str(), delta));
154 colSize.push_back(10);
155 colName.push_back("struct");
156 colInfo.push_back("Structure index.");
157 colSize.push_back(10);
158 colName.push_back("index");
159 colInfo.push_back("Property index.");
160 colSize.push_back(10);
161 colName.push_back("weight");
162 colInfo.push_back("Weight index.");
163 //colSize.push_back(10);
164 //colName.push_back("atom");
165 //colInfo.push_back("Atom index (starting with 1).");
166 //colSize.push_back(3);
167 //colName.push_back("xyz");
168 //colInfo.push_back("Force component (0 = x, 1 = y, 2 = z).");
169 colSize.push_back(24);
170 colName.push_back("analytic");
171 colInfo.push_back("Analytic weight derivatives.");
172 colSize.push_back(24);
173 colName.push_back("numeric");
174 colInfo.push_back("Numeric weight derivatives.");
175 appendLinesToFile(*(outFiles.at(k)),
176 createFileHeader(title,
177 colSize,
178 colName,
179 colInfo));
180 }
181 }
182
183 training.log << "\n";
184 training.log << " |";
185 for (auto k : pk) training.log << strpr(" %33s |", k.c_str());
186 training.log << "\n";
187 training.log << " |";
188 for (auto k : pk) training.log << " meanAbsError maxAbsError verdict |";
189 training.log << "\n";
190 training.log << "-----------------------";
191 for (auto k : pk) training.log << "------------------------------------";
192 training.log << "\n";
193
194
195 // Loop over all structures and compute symmetry functions and energies.
196 bool useForces = training.settingsKeywordExists("use_short_forces");
197 for (auto& s : training.structures)
198 {
199 s.calculateNeighborList(training.getMaxCutoffRadius());
200#ifdef N2P2_NO_SF_GROUPS
201 training.calculateSymmetryFunctions(s, useForces);
202#else
203 training.calculateSymmetryFunctionGroups(s, useForces);
204#endif
205 training.log << strpr("Configuration %6zu: ", s.index + 1);
206 outFileSummary << strpr("%10zu %10zu", s.index + 1, s.numAtoms);
207 for (auto k : pk)
208 {
209 vector<vector<double>> dPdc;
210 training.dPdc(k, s, dPdc);
211 vector<vector<double>> dPdcN;
212 training.dPdcN(k, s, dPdcN, delta);
213 meanAbsError[k] = 0.0;
214 maxAbsError[k] = 0.0;
215 size_t count = 0;
216 for (size_t i = 0; i < dPdc.size(); ++i)
217 {
218 for (size_t j = 0; j < dPdc.at(i).size(); ++j)
219 {
220 *(outFiles.at(k)) << strpr("%10zu %10zu %10zu "
221 "%24.16E %24.16E\n",
222 s.index + 1,
223 i + 1,
224 j + 1,
225 dPdc.at(i).at(j),
226 dPdcN.at(i).at(j));
227 outFiles.at(k)->flush();
228 double const error = fabs(dPdc.at(i).at(j) -
229 dPdcN.at(i).at(j));
230 meanAbsError.at(k) += error;
231 maxAbsError.at(k) = max(error, maxAbsError.at(k));
232 count++;
233 }
234 }
235 meanAbsError.at(k) /= count;
236 string swarn = "OK";
237 if (maxAbsError.at(k) > 100 * delta * delta)
238 {
239 swarn = "WARN";
240 warning++;
241 }
242 training.log << strpr(" %12.3E %12.3E %7s ",
243 meanAbsError.at(k),
244 maxAbsError.at(k),
245 swarn.c_str());
246 outFileSummary << strpr(" %16.8E %16.8E",
247 meanAbsError.at(k),
248 maxAbsError.at(k));
249 }
250 training.log << "\n";
251 outFileSummary << "\n";
252 s.reset();
253 }
254 if (normalize) training.toPhysicalUnits();
255
256 outFileSummary.close();
257 MPI_Barrier(MPI_COMM_WORLD);
258 if (myRank == 0) training.combineFiles("checkdw-summary.out");
259
260 for (auto k : pk)
261 {
262 outFiles.at(k)->close();
263 delete outFiles.at(k);
264 MPI_Barrier(MPI_COMM_WORLD);
265 if (myRank == 0) training.combineFiles(outFilesNames.at(k));
266 }
267
268 MPI_Allreduce(MPI_IN_PLACE, &warning, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
269 if (warning > 0)
270 {
271 training.log << "\n";
272 training.log << "IMPORTANT: Some warnings were issued. By default, this happens if the maximum\n"
273 " absolute error (\"maxAbsError\") is higher than 100 * delta².\n"
274 " However, this does NOT mean that analytic derivatives are incorrect!\n"
275 " Repeat this analysis with a different delta and check whether\n"
276 " the error scales with O(delta²). The prefactor for your system\n"
277 " could be higher than 10, hence, as long as there is a O(delta²)\n"
278 " scaling the analytic derivatives are probably correct.\n";
279 }
280
281 training.log << "*****************************************"
282 "**************************************\n";
283
284 myLog.close();
285
286 MPI_Finalize();
287
288 return 0;
289}
int distributeStructures(bool randomize, bool excludeRank0=false, std::string const &fileName="input.data")
Read data file and distribute structures among processors.
Definition: Dataset.cpp:724
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 toPhysicalUnits()
Switch all structures to physical units.
Definition: Dataset.cpp:963
void combineFiles(std::string filePrefix) const
Combine individual MPI proc files to one.
Definition: Dataset.cpp:1744
void toNormalizedUnits()
Switch all structures to normalized units.
Definition: Dataset.cpp:952
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
void initialize()
Write welcome message with version information.
Definition: Mode.cpp:55
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:212
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:1445
double getMaxCutoffRadius() const
Getter for Mode::maxCutoffRadius.
Definition: Mode.h:693
bool useNormalization() const
Check if normalization is enabled.
Definition: Mode.h:703
void calculateSymmetryFunctionGroups(Structure &structure, bool const derivatives)
Calculate all symmetry function groups for all atoms in given structure.
Definition: Mode.cpp:1561
virtual void setupSymmetryFunctionScaling(std::string const &fileName="scaling.data")
Set up symmetry function scaling from file.
Definition: Mode.cpp:712
bool settingsKeywordExists(std::string const &keyword) const
Check if keyword was found in settings file.
Definition: Mode.cpp:2193
Log log
Global log file.
Definition: Mode.h:593
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
Definition: Mode.cpp:1480
void setupSymmetryFunctionStatistics(bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
Set up symmetry function statistics collection.
Definition: Mode.cpp:1103
void loadSettingsFile(std::string const &fileName="input.nn")
Open settings file and load all keywords into memory.
Definition: Mode.cpp:161
Training methods.
Definition: Training.h:36
void dPdc(std::string property, Structure &structure, std::vector< std::vector< double > > &dEdc)
Compute derivatives of property with respect to weights.
Definition: Training.cpp:3307
void dPdcN(std::string property, Structure &structure, std::vector< std::vector< double > > &dEdc, double delta=1.0E-4)
Compute numeric derivatives of property with respect to weights.
Definition: Training.cpp:3379
std::vector< std::string > setupNumericDerivCheck()
Set up numeric weight derivatives check.
Definition: Training.cpp:1209
void setStage(std::size_t stage)
Set training stage (if multiple stages are needed for NNP type).
Definition: Training.cpp:379
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:111
void appendLinesToFile(ofstream &file, vector< string > const lines)
Append multiple lines of strings to open file stream.
Definition: utility.cpp:225

References nnp::appendLinesToFile(), nnp::Mode::calculateSymmetryFunctionGroups(), nnp::Mode::calculateSymmetryFunctions(), nnp::Dataset::combineFiles(), nnp::createFileHeader(), nnp::Dataset::distributeStructures(), nnp::Training::dPdc(), nnp::Training::dPdcN(), nnp::Mode::getMaxCutoffRadius(), nnp::Mode::initialize(), nnp::Mode::loadSettingsFile(), nnp::Mode::log, nnp::Log::registerStreamPointer(), nnp::Training::setStage(), nnp::Mode::settingsKeywordExists(), nnp::Mode::setupGeneric(), nnp::Dataset::setupMPI(), nnp::Mode::setupNeuralNetworkWeights(), nnp::Training::setupNumericDerivCheck(), nnp::Mode::setupSymmetryFunctionScaling(), nnp::Mode::setupSymmetryFunctionStatistics(), nnp::strpr(), nnp::Dataset::structures, nnp::Dataset::toNormalizedUnits(), nnp::Dataset::toPhysicalUnits(), nnp::Mode::useNormalization(), and nnp::Log::writeToStdout.

Here is the call graph for this function: