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
63 myLog.open((
strpr(
"nnp-norm2.log.%04d", myRank) + suffix).c_str());
72
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
83 {
84 throw runtime_error("ERROR: Normalization keywords found in settings, "
85 "please remove them first.\n");
86 }
87
89 {
90 throw runtime_error("ERROR: Normalization is only possible if forces "
91 "are used (keyword \"use_short_forces\").\n");
92 }
93
94
96
97
99
100
102
104 if (nnpType == Training::NNPType::SHORT_CHARGE_NN && stage == 2)
105 {
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
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).");
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";
158 {
159
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,
168#ifdef N2P2_NO_SF_GROUPS
170#else
172#endif
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";
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;
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
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
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");
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";
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.
void setupMPI()
Initialize MPI with MPI_COMM_WORLD.
std::vector< Structure > structures
All structures in this dataset.
void setupRandomNumberGenerator()
Initialize random number generator.
void combineFiles(std::string filePrefix) const
Combine individual MPI proc files to one.
void registerStreamPointer(std::ofstream *const &streamPointer)
Register new C++ ofstream pointer.
bool writeToStdout
Turn on/off output to stdout.
NNPType getNnpType() const
Getter for Mode::nnpType.
void initialize()
Write welcome message with version information.
double getMaxCutoffRadius() const
Getter for Mode::maxCutoffRadius.
void calculateAtomicNeuralNetworks(Structure &structure, bool const derivatives)
Calculate a single atomic neural network for a given atom and nn type.
void calculateEnergy(Structure &structure) const
Calculate potential energy for a given structure.
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.
double getEnergyWithOffset(Structure const &structure, bool ref=true) const
Add atomic energy offsets and return energy.
void setupGeneric(bool skipNormalize=false)
Combine multiple setup routines and provide a basic NNP setup.
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
void calculateForces(Structure &structure) const
Calculate forces 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 writeSettingsFile(std::ofstream *const &file) const
Write complete settings file.
void initializeWeights()
Initialize weights for all elements.
void writeWeights(std::string const &nnName, std::string const &fileNameFormat) const
Write weights to files (one file for each element).
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.