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
72 myLog.open(
strpr(
"nnp-checkdw.log.%04d", myRank).c_str());
86
88 training.
log <<
"*** ANALYTIC/NUMERIC WEIGHT DERIVATIVES C"
89 "HECK *********************************\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
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 }
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
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
164
165
166
167
168
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.");
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
198 {
200#ifdef N2P2_NO_SF_GROUPS
202#else
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 }
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.
void setupMPI()
Initialize MPI with MPI_COMM_WORLD.
std::vector< Structure > structures
All structures in this dataset.
void toPhysicalUnits()
Switch all structures to physical units.
void combineFiles(std::string filePrefix) const
Combine individual MPI proc files to one.
void toNormalizedUnits()
Switch all structures to normalized units.
void registerStreamPointer(std::ofstream *const &streamPointer)
Register new C++ ofstream pointer.
bool writeToStdout
Turn on/off output to stdout.
void initialize()
Write welcome message with version information.
void setupGeneric(std::string const &nnpDir="", bool skipNormalize=false, bool initialHardness=false)
Combine multiple setup routines and provide a basic NNP setup.
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.
double getMaxCutoffRadius() const
Getter for Mode::maxCutoffRadius.
bool useNormalization() const
Check if normalization is enabled.
void calculateSymmetryFunctionGroups(Structure &structure, bool const derivatives)
Calculate all symmetry function groups for all atoms in given structure.
virtual void setupSymmetryFunctionScaling(std::string const &fileName="scaling.data")
Set up symmetry function scaling from file.
bool settingsKeywordExists(std::string const &keyword) const
Check if keyword was found in settings file.
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
void setupSymmetryFunctionStatistics(bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
Set up symmetry function statistics collection.
void loadSettingsFile(std::string const &fileName="input.nn")
Open settings file and load all keywords into memory.
void dPdc(std::string property, Structure &structure, std::vector< std::vector< double > > &dEdc)
Compute derivatives of property with respect to weights.
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.
std::vector< std::string > setupNumericDerivCheck()
Set up numeric weight derivatives check.
void setStage(std::size_t stage)
Set training stage (if multiple stages are needed for NNP type).
string strpr(const char *format,...)
String version of printf function.
vector< string > createFileHeader(vector< string > const &title, vector< size_t > const &colSize, vector< string > const &colName, vector< string > const &colInfo, char const &commentChar)
void appendLinesToFile(ofstream &file, vector< string > const lines)
Append multiple lines of strings to open file stream.