n2p2 - A neural network potential package
Loading...
Searching...
No Matches
Mode.cpp
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#include "Mode.h"
18#include "NeuralNetwork.h"
19#include "utility.h"
20#include "version.h"
21#include <cmath>
22#ifdef _OPENMP
23#include <omp.h>
24#endif
25#include <Eigen/QR>
26#include <algorithm> // std::min, std::max, std::remove_if
27#include <cstdlib> // atoi, atof
28#include <fstream> // std::ifstream
29#ifndef N2P2_NO_SF_CACHE
30#include <map> // std::multimap
31#endif
32#include <limits> // std::numeric_limits
33#include <stdexcept> // std::runtime_error
34#include <utility> // std::piecewise_construct, std::forward_as_tuple
35
36using namespace std;
37using namespace nnp;
38using namespace Eigen;
39
41 normalize (false ),
43 numElements (0 ),
44 maxCutoffRadius (0.0 ),
45 cutoffAlpha (0.0 ),
46 meanEnergy (0.0 ),
47 convEnergy (1.0 ),
48 convLength (1.0 ),
49 convCharge (1.0 ),
50 fourPiEps (1.0 ),
51 ewaldSetup {},
52 kspaceGrid {}
53{
54}
55
57{
58 log << "\n";
59 log << "*****************************************"
60 "**************************************\n";
61 log << "\n";
62 log << "WELCOME TO n²p², A SOFTWARE PACKAGE FOR NEURAL NETWORK "
63 "POTENTIALS!\n";
64 log << "-------------------------------------------------------"
65 "-----------\n";
66 log << "\n";
67 log << "n²p² version (from git): " N2P2_GIT_VERSION "\n";
68 log << " (version.h): " N2P2_VERSION "\n";
69 log << "------------------------------------------------------------\n";
70 log << "Git branch : " N2P2_GIT_BRANCH "\n";
71 log << "Git revision : " N2P2_GIT_REV "\n";
72 log << "Compile date/time : " __DATE__ " " __TIME__ "\n";
73 log << "------------------------------------------------------------\n";
74 log << "\n";
75 log << "Features/Flags:\n";
76 log << "------------------------------------------------------------\n";
77 log << "Symmetry function groups : ";
78#ifndef N2P2_NO_SF_GROUPS
79 log << "enabled\n";
80#else
81 log << "disabled\n";
82#endif
83 log << "Symmetry function cache : ";
84#ifndef N2P2_NO_SF_CACHE
85 log << "enabled\n";
86#else
87 log << "disabled\n";
88#endif
89 log << "Timing function available : ";
90#ifndef N2P2_NO_TIME
91 log << "available\n";
92#else
93 log << "disabled\n";
94#endif
95 log << "Asymmetric polynomial SFs : ";
96#ifndef N2P2_NO_ASYM_POLY
97 log << "available\n";
98#else
99 log << "disabled\n";
100#endif
101 log << "SF low neighbor number check : ";
102#ifndef N2P2_NO_NEIGH_CHECK
103 log << "enabled\n";
104#else
105 log << "disabled\n";
106#endif
107 log << "SF derivative memory layout : ";
108#ifndef N2P2_FULL_SFD_MEMORY
109 log << "reduced\n";
110#else
111 log << "full\n";
112#endif
113 log << "MPI explicitly disabled : ";
114#ifndef N2P2_NO_MPI
115 log << "no\n";
116#else
117 log << "yes\n";
118#endif
119#ifdef _OPENMP
120 log << strpr("OpenMP threads : %d\n", omp_get_max_threads());
121#endif
122 log << "------------------------------------------------------------\n";
123 log << "\n";
124
125 log << "Please cite the following papers when publishing results "
126 "obtained with n²p²:\n";
127 log << "-----------------------------------------"
128 "--------------------------------------\n";
129 log << " * General citation for n²p² and the LAMMPS interface:\n";
130 log << "\n";
131 log << " Singraber, A.; Behler, J.; Dellago, C.\n";
132 log << " Library-Based LAMMPS Implementation of High-Dimensional\n";
133 log << " Neural Network Potentials.\n";
134 log << " J. Chem. Theory Comput. 2019 15 (3), 1827–1840.\n";
135 log << " https://doi.org/10.1021/acs.jctc.8b00770\n";
136 log << "-----------------------------------------"
137 "--------------------------------------\n";
138 log << " * Additionally, if you use the NNP training features of n²p²:\n";
139 log << "\n";
140 log << " Singraber, A.; Morawietz, T.; Behler, J.; Dellago, C.\n";
141 log << " Parallel Multistream Training of High-Dimensional Neural\n";
142 log << " Network Potentials.\n";
143 log << " J. Chem. Theory Comput. 2019, 15 (5), 3075–3092.\n";
144 log << " https://doi.org/10.1021/acs.jctc.8b01092\n";
145 log << "-----------------------------------------"
146 "--------------------------------------\n";
147 log << " * Additionally, if polynomial symmetry functions are used:\n";
148 log << "\n";
149 log << " Bircher, M. P.; Singraber, A.; Dellago, C.\n";
150 log << " Improved Description of Atomic Environments Using Low-Cost\n";
151 log << " Polynomial Functions with Compact Support.\n";
152 log << " arXiv:2010.14414 [cond-mat, physics:physics] 2020.\n";
153 log << " https://arxiv.org/abs/2010.14414\n";
154
155
156 log << "*****************************************"
157 "**************************************\n";
158
159 return;
160}
161
162void Mode::loadSettingsFile(string const& fileName)
163{
164 log << "\n";
165 log << "*** SETUP: SETTINGS FILE ****************"
166 "**************************************\n";
167 log << "\n";
168
169 size_t numCriticalProblems = settings.loadFile(fileName);
170 log << settings.info();
171 if (numCriticalProblems > 0)
172 {
173 throw runtime_error(strpr("ERROR: %zu critical problem(s) were found "
174 "in settings file.\n", numCriticalProblems));
175 }
176
177 if (settings.keywordExists("nnp_type"))
178 {
179 string nnpTypeString = settings["nnp_type"];
180 if (nnpTypeString == "2G-HDNNP") nnpType = NNPType::HDNNP_2G;
181 else if (nnpTypeString == "4G-HDNNP") nnpType = NNPType::HDNNP_4G;
182 else if (nnpTypeString == "Q-HDNNP") nnpType = NNPType::HDNNP_Q;
183 else nnpType = (NNPType)atoi(settings["nnp_type"].c_str());
184 }
185
187 {
188 log << "This settings file defines a short-range NNP (2G-HDNNP).\n";
189 }
190 else if (nnpType == NNPType::HDNNP_4G)
191 {
192 log << "This settings file defines a NNP with electrostatics and\n"
193 "non-local charge transfer (4G-HDNNP).\n";
194 }
195 else if (nnpType == NNPType::HDNNP_Q)
196 {
197 log << "This settings file defines a short-range NNP similar to\n"
198 "4G-HDNNP with additional charge NN but with neither\n"
199 "electrostatics nor global charge equilibration\n"
200 "(method by M. Bircher).\n";
201 }
202 else
203 {
204 throw runtime_error("ERROR: Unknown NNP type.\n");
205 }
206
207 log << "*****************************************"
208 "**************************************\n";
209
210 return;
211}
212
213void Mode::setupGeneric(string const& nnpDir,
214 bool skipNormalize,
215 bool initialHardness)
216{
217 if (!skipNormalize) setupNormalization();
221 setupElectrostatics(initialHardness, nnpDir);
222 setupCutoff();
225#ifndef N2P2_FULL_SFD_MEMORY
227#endif
228#ifndef N2P2_NO_SF_CACHE
230#endif
231#ifndef N2P2_NO_SF_GROUPS
233#endif
235
236 return;
237}
238
239void Mode::setupNormalization(bool standalone)
240{
241 if (standalone)
242 {
243 log << "\n";
244 log << "*** SETUP: NORMALIZATION ****************"
245 "**************************************\n";
246 log << "\n";
247 }
248
249 if (settings.keywordExists("mean_energy") &&
250 settings.keywordExists("conv_energy") &&
251 settings.keywordExists("conv_length"))
252 {
253 normalize = true;
254 meanEnergy = atof(settings["mean_energy"].c_str());
255 convEnergy = atof(settings["conv_energy"].c_str());
256 convLength = atof(settings["conv_length"].c_str());
257
258 log << "Data set normalization is used.\n";
259 log << strpr("Mean energy per atom : %24.16E\n", meanEnergy);
260 log << strpr("Conversion factor energy : %24.16E\n", convEnergy);
261 log << strpr("Conversion factor length : %24.16E\n", convLength);
262 if (settings.keywordExists("conv_charge"))
263 {
264 convCharge = atof(settings["conv_charge"].c_str());
265 log << strpr("Conversion factor charge : %24.16E\n",
266 convCharge);
267 }
268 if (settings.keywordExists("atom_energy"))
269 {
270 log << "\n";
271 log << "Atomic energy offsets are used in addition to"
272 " data set normalization.\n";
273 log << "Offsets will be subtracted from reference energies BEFORE"
274 " normalization is applied.\n";
275 }
276 }
277 else if ((!settings.keywordExists("mean_energy")) &&
278 (!settings.keywordExists("conv_energy")) &&
279 (!settings.keywordExists("conv_length")) &&
280 (!settings.keywordExists("conv_charge")))
281 {
282 normalize = false;
283 log << "Data set normalization is not used.\n";
284 }
285 else
286 {
287 throw runtime_error("ERROR: Incorrect usage of normalization"
288 " keywords.\n"
289 " Use all or none of \"mean_energy\", \"conv_energy\""
290 " and \"conv_length\".\n");
291 }
292
293 if (standalone)
294 {
295 log << "*****************************************"
296 "**************************************\n";
297 }
298
299 return;
300}
301
303{
304 log << "\n";
305 log << "*** SETUP: ELEMENT MAP ******************"
306 "**************************************\n";
307 log << "\n";
308
309 elementMap.registerElements(settings["elements"]);
310 log << strpr("Number of element strings found: %d\n", elementMap.size());
311 for (size_t i = 0; i < elementMap.size(); ++i)
312 {
313 log << strpr("Element %2zu: %2s (%3zu)\n", i, elementMap[i].c_str(),
314 elementMap.atomicNumber(i));
315 }
316
317 log << "*****************************************"
318 "**************************************\n";
319
320 return;
321}
322
324{
325 log << "\n";
326 log << "*** SETUP: ELEMENTS *********************"
327 "**************************************\n";
328 log << "\n";
329
330 numElements = (size_t)atoi(settings["number_of_elements"].c_str());
331 if (numElements != elementMap.size())
332 {
333 throw runtime_error("ERROR: Inconsistent number of elements.\n");
334 }
335 log << strpr("Number of elements is consistent: %zu\n", numElements);
336
337 for (size_t i = 0; i < numElements; ++i)
338 {
339 elements.push_back(Element(i, elementMap));
340 }
341
342 if (settings.keywordExists("atom_energy"))
343 {
344 settings::Settings::KeyRange r = settings.getValues("atom_energy");
345 for (settings::Settings::KeyMap::const_iterator it = r.first;
346 it != r.second; ++it)
347 {
348 vector<string> args = split(reduce(it->second.first));
349 size_t element = elementMap[args.at(0)];
350 elements.at(element).
351 setAtomicEnergyOffset(atof(args.at(1).c_str()));
352 }
353 }
354 log << "Atomic energy offsets per element:\n";
355 for (size_t i = 0; i < elementMap.size(); ++i)
356 {
357 log << strpr("Element %2zu: %16.8E\n",
358 i, elements.at(i).getAtomicEnergyOffset());
359 }
360 log << "Energy offsets are automatically subtracted from reference "
361 "energies.\n";
362
363 log << "*****************************************"
364 "**************************************\n";
365
366 return;
367}
368
369void Mode::setupElectrostatics(bool initialHardness,
370 string directoryPrefix,
371 string fileNameFormat)
372{
373 log << "\n";
374 log << "*** SETUP: ELECTROSTATICS ***************"
375 "**************************************\n";
376 log << "\n";
377
378 // Atomic hardness.
379 if (initialHardness)
380 {
381 settings::Settings::KeyRange r = settings.getValues("initial_hardness");
382 for (settings::Settings::KeyMap::const_iterator it = r.first;
383 it != r.second; ++it)
384 {
385 vector<string> args = split(reduce(it->second.first));
386 size_t element = elementMap[args.at(0)];
387 double hardness = atof(args.at(1).c_str());
388 if (normalize) hardness = normalized("hardness", hardness);
389 elements.at(element).setHardness(hardness);
390 }
391 for (size_t i = 0; i < numElements; ++i)
392 {
393 double hardness = elements.at(i).getHardness();
394 if (normalize) hardness = physical("hardness", hardness);
395 log << strpr("Initial atomic hardness for element %2s: %16.8E\n",
396 elements.at(i).getSymbol().c_str(),
397 hardness);
398 }
399 }
400 else
401 {
402 string actualFileNameFormat = directoryPrefix + fileNameFormat;
403 log << strpr("Atomic hardness file name format: %s\n",
404 actualFileNameFormat.c_str());
405 for (size_t i = 0; i < numElements; ++i)
406 {
407 string fileName = strpr(actualFileNameFormat.c_str(),
408 elements.at(i).getAtomicNumber());
409 log << strpr("Atomic hardness for element %2s from file %s: ",
410 elements.at(i).getSymbol().c_str(),
411 fileName.c_str());
412 vector<double> const data = readColumnsFromFile(fileName,
413 {0}).at(0);
414 if (data.size() != 1)
415 {
416 throw runtime_error("ERROR: Atomic hardness data is "
417 "inconsistent.\n");
418 }
419 double hardness = data.at(0);
420 if (normalize) hardness = normalized("hardness", hardness);
421 elements.at(i).setHardness(hardness);
422 if (normalize) hardness = physical("hardness", hardness);
423 log << strpr("%16.8E\n", hardness);
424 }
425 }
426 log << "\n";
427
428 // Gaussian width of charges.
429 double maxQsigma = 0.0;
430 settings::Settings::KeyRange r = settings.getValues("fixed_gausswidth");
431 for (settings::Settings::KeyMap::const_iterator it = r.first;
432 it != r.second; ++it)
433 {
434 vector<string> args = split(reduce(it->second.first));
435 size_t element = elementMap[args.at(0)];
436 double qsigma = atof(args.at(1).c_str());
437 if (normalize) qsigma *= convLength;
438 elements.at(element).setQsigma(qsigma);
439
440 maxQsigma = max(qsigma, maxQsigma);
441 }
442 log << "Gaussian width of charge distribution per element:\n";
443 for (size_t i = 0; i < elementMap.size(); ++i)
444 {
445 double qsigma = elements.at(i).getQsigma();
446 if (normalize) qsigma /= convLength;
447 log << strpr("Element %2zu: %16.8E\n",
448 i, qsigma);
449 }
450 log << "\n";
451
452 // K-Space Solver
453 if (settings.keywordExists("kspace_solver"))
454 {
455 string kspace_solver_string = settings["kspace_solver"];
456 if (kspace_solver_string == "ewald")
457 {
458 kspaceGrid.kspaceSolver = KSPACESolver::EWALD_SUM;
459 }
460 else if (kspace_solver_string == "pppm")
461 {
462 kspaceGrid.kspaceSolver = KSPACESolver::PPPM;
463 }
464 else if (kspace_solver_string == "ewald_lammps")
465 {
467 }
468 }
469 else kspaceGrid.kspaceSolver = KSPACESolver::EWALD_SUM;
470
471 // Ewald truncation error method.
472 if (settings.keywordExists("ewald_truncation_error_method"))
473 {
474 string truncation_method_string =
475 settings["ewald_truncation_error_method"];
476 if (truncation_method_string == "0")
477 {
479 }
480 else if (truncation_method_string == "1")
481 {
483 }
484 }
485 else ewaldSetup.setTruncMethod(EWALDTruncMethod::JACKSON_CATLOW);
486
487 // Ewald precision.
488 if (settings.keywordExists("ewald_prec"))
489 {
490 vector<string> args = split(settings["ewald_prec"]);
491 ewaldSetup.readFromArgs(args);
492 ewaldSetup.setMaxQSigma(maxQsigma);
493 if (normalize) ewaldSetup.toNormalizedUnits(convEnergy, convLength);
494 }
495 else if (ewaldSetup.getTruncMethod() == EWALDTruncMethod::KOLAFA_PERRAM)
496 {
497 throw runtime_error("ERROR: ewald_truncation_error_method 1 requires "
498 "ewald_prec.");
499 }
500 log << strpr("Ewald truncation error method: %16d\n",
501 ewaldSetup.getTruncMethod());
502 log << strpr("Ewald precision: %16.8E\n",
503 ewaldSetup.getPrecision());
504 if (ewaldSetup.getTruncMethod() == EWALDTruncMethod::KOLAFA_PERRAM)
505 {
506 log << strpr("Ewald expected maximum charge: %16.8E\n",
507 ewaldSetup.getMaxCharge());
508 }
509 log << "\n";
510
511 // 4 * pi * epsilon
512 if (settings.keywordExists("four_pi_epsilon"))
513 fourPiEps = atof(settings["four_pi_epsilon"].c_str());
514 if (normalize)
515 fourPiEps *= pow(convCharge, 2) / (convLength * convEnergy);
516
517 // Screening function.
518 if (settings.keywordExists("screen_electrostatics"))
519 {
520 vector<string> args = split(settings["screen_electrostatics"]);
521 double inner = atof(args.at(0).c_str());
522 double outer = atof(args.at(1).c_str());
523 string type = "c"; // Default core function is cosine.
524 if (args.size() > 2) type = args.at(2);
525 screeningFunction.setInnerOuter(inner, outer);
526 screeningFunction.setCoreFunction(type);
527 for (auto s : screeningFunction.info()) log << s;
528 if (normalize) screeningFunction.changeLengthUnits(convLength);
529 }
530 else
531 {
532 // Set screening function in such way that it will always be 1.0.
533 // This may not be very efficient, may optimize later.
534 screeningFunction.setInnerOuter(-2.0, -1.0);
535 log << "Screening function not used.\n";
536 }
537
538 log << "*****************************************"
539 "**************************************\n";
540
541 return;
542}
543
545{
546 log << "\n";
547 log << "*** SETUP: CUTOFF FUNCTIONS *************"
548 "**************************************\n";
549 log << "\n";
550
551 vector<string> args = split(settings["cutoff_type"]);
552
553 cutoffType = (CutoffFunction::CutoffType) atoi(args.at(0).c_str());
554 if (args.size() > 1)
555 {
556 cutoffAlpha = atof(args.at(1).c_str());
557 if (0.0 < cutoffAlpha && cutoffAlpha >= 1.0)
558 {
559 throw invalid_argument("ERROR: 0 <= alpha < 1.0 is required.\n");
560 }
561 }
562 log << strpr("Parameter alpha for inner cutoff: %f\n", cutoffAlpha);
563 log << "Inner cutoff = Symmetry function cutoff * alpha\n";
564
565 log << "Equal cutoff function type for all symmetry functions:\n";
567 {
568 log << strpr("CutoffFunction::CT_COS (%d)\n", cutoffType);
569 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
570 log << "f(x) = 1/2 * (cos(pi*x) + 1)\n";
571 }
573 {
574 log << strpr("CutoffFunction::CT_TANHU (%d)\n", cutoffType);
575 log << "f(r) = tanh^3(1 - r/rc)\n";
576 if (cutoffAlpha > 0.0)
577 {
578 log << "WARNING: Inner cutoff parameter not used in combination"
579 " with this cutoff function.\n";
580 }
581 }
583 {
584 log << strpr("CutoffFunction::CT_TANH (%d)\n", cutoffType);
585 log << "f(r) = c * tanh^3(1 - r/rc), f(0) = 1\n";
586 if (cutoffAlpha > 0.0)
587 {
588 log << "WARNING: Inner cutoff parameter not used in combination"
589 " with this cutoff function.\n";
590 }
591 }
593 {
594 log << strpr("CutoffFunction::CT_POLY1 (%d)\n", cutoffType);
595 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
596 log << "f(x) = (2x - 3)x^2 + 1\n";
597 }
599 {
600 log << strpr("CutoffFunction::CT_POLY2 (%d)\n", cutoffType);
601 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
602 log << "f(x) = ((15 - 6x)x - 10)x^3 + 1\n";
603 }
605 {
606 log << strpr("CutoffFunction::CT_POLY3 (%d)\n", cutoffType);
607 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
608 log << "f(x) = (x(x(20x - 70) + 84) - 35)x^4 + 1\n";
609 }
611 {
612 log << strpr("CutoffFunction::CT_POLY4 (%d)\n", cutoffType);
613 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
614 log << "f(x) = (x(x((315 - 70x)x - 540) + 420) - 126)x^5 + 1\n";
615 }
617 {
618 log << strpr("CutoffFunction::CT_EXP (%d)\n", cutoffType);
619 log << "x := (r - rc * alpha) / (rc - rc * alpha)\n";
620 log << "f(x) = exp(-1 / 1 - x^2)\n";
621 }
623 {
624 log << strpr("CutoffFunction::CT_HARD (%d)\n", cutoffType);
625 log << "f(r) = 1\n";
626 log << "WARNING: Hard cutoff used!\n";
627 }
628 else
629 {
630 throw invalid_argument("ERROR: Unknown cutoff type.\n");
631 }
632
633 log << "*****************************************"
634 "**************************************\n";
635
636 return;
637}
638
640{
641 log << "\n";
642 log << "*** SETUP: SYMMETRY FUNCTIONS ***********"
643 "**************************************\n";
644 log << "\n";
645
646 settings::Settings::KeyRange r = settings.getValues("symfunction_short");
647 for (settings::Settings::KeyMap::const_iterator it = r.first; it != r.second; ++it)
648 {
649 vector<string> args = split(reduce(it->second.first));
650 size_t element = elementMap[args.at(0)];
651
652 elements.at(element).addSymmetryFunction(it->second.first,
653 it->second.second);
654 }
655
656 log << "Abbreviations:\n";
657 log << "--------------\n";
658 log << "ind .... Symmetry function index.\n";
659 log << "ec ..... Central atom element.\n";
660 log << "tp ..... Symmetry function type.\n";
661 log << "sbtp ... Symmetry function subtype (e.g. cutoff type).\n";
662 log << "e1 ..... Neighbor 1 element.\n";
663 log << "e2 ..... Neighbor 2 element.\n";
664 log << "eta .... Gaussian width eta.\n";
665 log << "rs/rl... Shift distance of Gaussian or left cutoff radius "
666 "for polynomial.\n";
667 log << "angl.... Left cutoff angle for polynomial.\n";
668 log << "angr.... Right cutoff angle for polynomial.\n";
669 log << "la ..... Angle prefactor lambda.\n";
670 log << "zeta ... Angle term exponent zeta.\n";
671 log << "rc ..... Cutoff radius / right cutoff radius for polynomial.\n";
672 log << "a ...... Free parameter alpha (e.g. cutoff alpha).\n";
673 log << "ln ..... Line number in settings file.\n";
674 log << "\n";
675 maxCutoffRadius = 0.0;
676 for (vector<Element>::iterator it = elements.begin();
677 it != elements.end(); ++it)
678 {
679 if (normalize) it->changeLengthUnitSymmetryFunctions(convLength);
680 it->sortSymmetryFunctions();
681 maxCutoffRadius = max(it->getMaxCutoffRadius(), maxCutoffRadius);
682 it->setCutoffFunction(cutoffType, cutoffAlpha);
683 log << strpr("Short range atomic symmetry functions element %2s :\n",
684 it->getSymbol().c_str());
685 log << "--------------------------------------------------"
686 "-----------------------------------------------\n";
687 log << " ind ec tp sbtp e1 e2 eta rs/rl "
688 "rc angl angr la zeta a ln\n";
689 log << "--------------------------------------------------"
690 "-----------------------------------------------\n";
691 log << it->infoSymmetryFunctionParameters();
692 log << "--------------------------------------------------"
693 "-----------------------------------------------\n";
694 }
695 minNeighbors.clear();
696 minNeighbors.resize(numElements, 0);
697 minCutoffRadius.clear();
699 for (size_t i = 0; i < numElements; ++i)
700 {
701 minNeighbors.at(i) = elements.at(i).getMinNeighbors();
702 minCutoffRadius.at(i) = elements.at(i).getMinCutoffRadius();
703 log << strpr("Minimum cutoff radius for element %2s: %f\n",
704 elements.at(i).getSymbol().c_str(),
706 }
707 log << strpr("Maximum cutoff radius (global) : %f\n",
709
710 log << "*****************************************"
711 "**************************************\n";
712
713 return;
714}
715
717{
718 log << "\n";
719 log << "*** SETUP: SYMMETRY FUNCTION SCALING ****"
720 "**************************************\n";
721 log << "\n";
722
723 log << "No scaling for symmetry functions.\n";
724 for (vector<Element>::iterator it = elements.begin();
725 it != elements.end(); ++it)
726 {
727 it->setScalingNone();
728 }
729
730 log << "*****************************************"
731 "**************************************\n";
732
733 return;
734}
735
736void Mode::setupSymmetryFunctionScaling(string const& fileName)
737{
738 log << "\n";
739 log << "*** SETUP: SYMMETRY FUNCTION SCALING ****"
740 "**************************************\n";
741 log << "\n";
742
743 log << "Equal scaling type for all symmetry functions:\n";
744 if ( ( settings.keywordExists("scale_symmetry_functions" ))
745 && (!settings.keywordExists("center_symmetry_functions")))
746 {
748 log << strpr("Scaling type::ST_SCALE (%d)\n", scalingType);
749 log << "Gs = Smin + (Smax - Smin) * (G - Gmin) / (Gmax - Gmin)\n";
750 }
751 else if ( (!settings.keywordExists("scale_symmetry_functions" ))
752 && ( settings.keywordExists("center_symmetry_functions")))
753 {
755 log << strpr("Scaling type::ST_CENTER (%d)\n", scalingType);
756 log << "Gs = G - Gmean\n";
757 }
758 else if ( ( settings.keywordExists("scale_symmetry_functions" ))
759 && ( settings.keywordExists("center_symmetry_functions")))
760 {
762 log << strpr("Scaling type::ST_SCALECENTER (%d)\n", scalingType);
763 log << "Gs = Smin + (Smax - Smin) * (G - Gmean) / (Gmax - Gmin)\n";
764 }
765 else if (settings.keywordExists("scale_symmetry_functions_sigma"))
766 {
768 log << strpr("Scaling type::ST_SCALESIGMA (%d)\n", scalingType);
769 log << "Gs = Smin + (Smax - Smin) * (G - Gmean) / Gsigma\n";
770 }
771 else
772 {
774 log << strpr("Scaling type::ST_NONE (%d)\n", scalingType);
775 log << "Gs = G\n";
776 log << "WARNING: No symmetry function scaling!\n";
777 }
778
779 double Smin = 0.0;
780 double Smax = 0.0;
784 {
785 if (settings.keywordExists("scale_min_short"))
786 {
787 Smin = atof(settings["scale_min_short"].c_str());
788 }
789 else
790 {
791 log << "WARNING: Keyword \"scale_min_short\" not found.\n";
792 log << " Default value for Smin = 0.0.\n";
793 Smin = 0.0;
794 }
795
796 if (settings.keywordExists("scale_max_short"))
797 {
798 Smax = atof(settings["scale_max_short"].c_str());
799 }
800 else
801 {
802 log << "WARNING: Keyword \"scale_max_short\" not found.\n";
803 log << " Default value for Smax = 1.0.\n";
804 Smax = 1.0;
805 }
806
807 log << strpr("Smin = %f\n", Smin);
808 log << strpr("Smax = %f\n", Smax);
809 }
810
811 log << strpr("Symmetry function scaling statistics from file: %s\n",
812 fileName.c_str());
813 log << "-----------------------------------------"
814 "--------------------------------------\n";
815 ifstream file;
816 file.open(fileName.c_str());
817 if (!file.is_open())
818 {
819 throw runtime_error("ERROR: Could not open file: \"" + fileName
820 + "\".\n");
821 }
822 string line;
823 vector<string> lines;
824 while (getline(file, line))
825 {
826 if (line.at(0) != '#') lines.push_back(line);
827 }
828 file.close();
829
830 log << "\n";
831 log << "Abbreviations:\n";
832 log << "--------------\n";
833 log << "ind ..... Symmetry function index.\n";
834 log << "min ..... Minimum symmetry function value.\n";
835 log << "max ..... Maximum symmetry function value.\n";
836 log << "mean .... Mean symmetry function value.\n";
837 log << "sigma ... Standard deviation of symmetry function values.\n";
838 log << "sf ...... Scaling factor for derivatives.\n";
839 log << "Smin .... Desired minimum scaled symmetry function value.\n";
840 log << "Smax .... Desired maximum scaled symmetry function value.\n";
841 log << "t ....... Scaling type.\n";
842 log << "\n";
843 for (vector<Element>::iterator it = elements.begin();
844 it != elements.end(); ++it)
845 {
846 it->setScaling(scalingType, lines, Smin, Smax);
847 log << strpr("Scaling data for symmetry functions element %2s :\n",
848 it->getSymbol().c_str());
849 log << "-----------------------------------------"
850 "--------------------------------------\n";
851 log << " ind min max mean sigma sf Smin Smax t\n";
852 log << "-----------------------------------------"
853 "--------------------------------------\n";
854 log << it->infoSymmetryFunctionScaling();
855 log << "-----------------------------------------"
856 "--------------------------------------\n";
857 lines.erase(lines.begin(), lines.begin() + it->numSymmetryFunctions());
858 }
859
860 log << "*****************************************"
861 "**************************************\n";
862
863 return;
864}
865
867{
868 log << "\n";
869 log << "*** SETUP: SYMMETRY FUNCTION GROUPS *****"
870 "**************************************\n";
871 log << "\n";
872
873 log << "Abbreviations:\n";
874 log << "--------------\n";
875 log << "ind .... Symmetry function index.\n";
876 log << "ec ..... Central atom element.\n";
877 log << "tp ..... Symmetry function type.\n";
878 log << "sbtp ... Symmetry function subtype (e.g. cutoff type).\n";
879 log << "e1 ..... Neighbor 1 element.\n";
880 log << "e2 ..... Neighbor 2 element.\n";
881 log << "eta .... Gaussian width eta.\n";
882 log << "rs/rl... Shift distance of Gaussian or left cutoff radius "
883 "for polynomial.\n";
884 log << "angl.... Left cutoff angle for polynomial.\n";
885 log << "angr.... Right cutoff angle for polynomial.\n";
886 log << "la ..... Angle prefactor lambda.\n";
887 log << "zeta ... Angle term exponent zeta.\n";
888 log << "rc ..... Cutoff radius / right cutoff radius for polynomial.\n";
889 log << "a ...... Free parameter alpha (e.g. cutoff alpha).\n";
890 log << "ln ..... Line number in settings file.\n";
891 log << "mi ..... Member index.\n";
892 log << "sfi .... Symmetry function index.\n";
893 log << "e ...... Recalculate exponential term.\n";
894 log << "\n";
895 for (vector<Element>::iterator it = elements.begin();
896 it != elements.end(); ++it)
897 {
898 it->setupSymmetryFunctionGroups();
899 log << strpr("Short range atomic symmetry function groups "
900 "element %2s :\n", it->getSymbol().c_str());
901 log << "------------------------------------------------------"
902 "----------------------------------------------------\n";
903 log << " ind ec tp sbtp e1 e2 eta rs/rl "
904 "rc angl angr la zeta a ln mi sfi e\n";
905 log << "------------------------------------------------------"
906 "----------------------------------------------------\n";
907 log << it->infoSymmetryFunctionGroups();
908 log << "------------------------------------------------------"
909 "----------------------------------------------------\n";
910 }
911
912 log << "*****************************************"
913 "**************************************\n";
914
915 return;
916}
917
919{
920 log << "\n";
921 log << "*** SETUP: SYMMETRY FUNCTION MEMORY *****"
922 "**************************************\n";
923 log << "\n";
924
925 for (auto& e : elements)
926 {
927 e.setupSymmetryFunctionMemory();
928 vector<size_t> symmetryFunctionNumTable
929 = e.getSymmetryFunctionNumTable();
930 vector<vector<size_t>> symmetryFunctionTable
931 = e.getSymmetryFunctionTable();
932 log << strpr("Symmetry function derivatives memory table "
933 "for element %2s :\n", e.getSymbol().c_str());
934 log << "-----------------------------------------"
935 "--------------------------------------\n";
936 log << "Relevant symmetry functions for neighbors with element:\n";
937 for (size_t i = 0; i < numElements; ++i)
938 {
939 log << strpr("- %2s: %4zu of %4zu (%5.1f %)\n",
940 elementMap[i].c_str(),
941 symmetryFunctionNumTable.at(i),
942 e.numSymmetryFunctions(),
943 (100.0 * symmetryFunctionNumTable.at(i))
944 / e.numSymmetryFunctions());
945 if (verbose)
946 {
947 log << "-----------------------------------------"
948 "--------------------------------------\n";
949 for (auto isf : symmetryFunctionTable.at(i))
950 {
951 SymFnc const& sf = e.getSymmetryFunction(isf);
952 log << sf.parameterLine();
953 }
954 log << "-----------------------------------------"
955 "--------------------------------------\n";
956 }
957 }
958 log << "-----------------------------------------"
959 "--------------------------------------\n";
960 }
961 if (verbose)
962 {
963 for (auto& e : elements)
964 {
965 log << strpr("%2s - symmetry function per-element index table:\n",
966 e.getSymbol().c_str());
967 log << "-----------------------------------------"
968 "--------------------------------------\n";
969 log << " ind";
970 for (size_t i = 0; i < numElements; ++i)
971 {
972 log << strpr(" %4s", elementMap[i].c_str());
973 }
974 log << "\n";
975 log << "-----------------------------------------"
976 "--------------------------------------\n";
977 for (size_t i = 0; i < e.numSymmetryFunctions(); ++i)
978 {
979 SymFnc const& sf = e.getSymmetryFunction(i);
980 log << strpr("%4zu", sf.getIndex() + 1);
981 vector<size_t> indexPerElement = sf.getIndexPerElement();
982 for (auto ipe : sf.getIndexPerElement())
983 {
984 if (ipe == numeric_limits<size_t>::max())
985 {
986 log << strpr(" ");
987 }
988 else
989 {
990 log << strpr(" %4zu", ipe + 1);
991 }
992 }
993 log << "\n";
994 }
995 log << "-----------------------------------------"
996 "--------------------------------------\n";
997 }
998
999 }
1000
1001 log << "*****************************************"
1002 "**************************************\n";
1003
1004 return;
1005}
1006
1007#ifndef N2P2_NO_SF_CACHE
1009{
1010 log << "\n";
1011 log << "*** SETUP: SYMMETRY FUNCTION CACHE ******"
1012 "**************************************\n";
1013 log << "\n";
1014
1015 for (size_t i = 0; i < numElements; ++i)
1016 {
1017 using SFCacheList = Element::SFCacheList;
1018 vector<vector<SFCacheList>> cacheLists(numElements);
1019 Element& e = elements.at(i);
1020 for (size_t j = 0; j < e.numSymmetryFunctions(); ++j)
1021 {
1022 SymFnc const& s = e.getSymmetryFunction(j);
1023 for (auto identifier : s.getCacheIdentifiers())
1024 {
1025 size_t ne = atoi(split(identifier)[0].c_str());
1026 bool unknown = true;
1027 for (auto& c : cacheLists.at(ne))
1028 {
1029 if (identifier == c.identifier)
1030 {
1031 c.indices.push_back(s.getIndex());
1032 unknown = false;
1033 break;
1034 }
1035 }
1036 if (unknown)
1037 {
1038 cacheLists.at(ne).push_back(SFCacheList());
1039 cacheLists.at(ne).back().element = ne;
1040 cacheLists.at(ne).back().identifier = identifier;
1041 cacheLists.at(ne).back().indices.push_back(s.getIndex());
1042 }
1043 }
1044 }
1045 if (verbose)
1046 {
1047 log << strpr("Multiple cache identifiers for element %2s:\n\n",
1048 e.getSymbol().c_str());
1049 }
1050 double cacheUsageMean = 0.0;
1051 size_t cacheCount = 0;
1052 for (size_t j = 0; j < numElements; ++j)
1053 {
1054 if (verbose)
1055 {
1056 log << strpr("Neighbor %2s:\n", elementMap[j].c_str());
1057 }
1058 vector<SFCacheList>& c = cacheLists.at(j);
1059 c.erase(remove_if(c.begin(),
1060 c.end(),
1061 [](SFCacheList l)
1062 {
1063 return l.indices.size() <= 1;
1064 }), c.end());
1065 cacheCount += c.size();
1066 for (size_t k = 0; k < c.size(); ++k)
1067 {
1068 cacheUsageMean += c.at(k).indices.size();
1069 if (verbose)
1070 {
1071 log << strpr("Cache %zu, Identifier \"%s\", "
1072 "Symmetry functions",
1073 k, c.at(k).identifier.c_str());
1074 for (auto si : c.at(k).indices)
1075 {
1076 log << strpr(" %zu", si);
1077 }
1078 log << "\n";
1079 }
1080 }
1081 }
1082 e.setCacheIndices(cacheLists);
1083 //for (size_t j = 0; j < e.numSymmetryFunctions(); ++j)
1084 //{
1085 // SymFnc const& sf = e.getSymmetryFunction(j);
1086 // auto indices = sf.getCacheIndices();
1087 // size_t count = 0;
1088 // for (size_t k = 0; k < numElements; ++k)
1089 // {
1090 // count += indices.at(k).size();
1091 // }
1092 // if (count > 0)
1093 // {
1094 // log << strpr("SF %4zu:\n", sf.getIndex());
1095 // }
1096 // for (size_t k = 0; k < numElements; ++k)
1097 // {
1098 // if (indices.at(k).size() > 0)
1099 // {
1100 // log << strpr("- Neighbor %2s:", elementMap[k].c_str());
1101 // for (size_t l = 0; l < indices.at(k).size(); ++l)
1102 // {
1103 // log << strpr(" %zu", indices.at(k).at(l));
1104 // }
1105 // log << "\n";
1106 // }
1107 // }
1108 //}
1109 cacheUsageMean /= cacheCount;
1110 log << strpr("Element %2s: in total %zu caches, "
1111 "used %3.2f times on average.\n",
1112 e.getSymbol().c_str(), cacheCount, cacheUsageMean);
1113 if (verbose)
1114 {
1115 log << "-----------------------------------------"
1116 "--------------------------------------\n";
1117 }
1118 }
1119
1120 log << "*****************************************"
1121 "**************************************\n";
1122
1123 return;
1124}
1125#endif
1126
1127void Mode::setupSymmetryFunctionStatistics(bool collectStatistics,
1128 bool collectExtrapolationWarnings,
1129 bool writeExtrapolationWarnings,
1130 bool stopOnExtrapolationWarnings)
1131{
1132 log << "\n";
1133 log << "*** SETUP: SYMMETRY FUNCTION STATISTICS *"
1134 "**************************************\n";
1135 log << "\n";
1136
1137 log << "Equal symmetry function statistics for all elements.\n";
1138 log << strpr("Collect min/max/mean/sigma : %d\n",
1139 (int)collectStatistics);
1140 log << strpr("Collect extrapolation warnings : %d\n",
1141 (int)collectExtrapolationWarnings);
1142 log << strpr("Write extrapolation warnings immediately to stderr: %d\n",
1143 (int)writeExtrapolationWarnings);
1144 log << strpr("Halt on any extrapolation warning : %d\n",
1145 (int)stopOnExtrapolationWarnings);
1146 for (vector<Element>::iterator it = elements.begin();
1147 it != elements.end(); ++it)
1148 {
1149 it->statistics.collectStatistics = collectStatistics;
1150 it->statistics.collectExtrapolationWarnings =
1151 collectExtrapolationWarnings;
1152 it->statistics.writeExtrapolationWarnings = writeExtrapolationWarnings;
1153 it->statistics.stopOnExtrapolationWarnings =
1154 stopOnExtrapolationWarnings;
1155 }
1156
1157 checkExtrapolationWarnings = collectStatistics
1158 || collectExtrapolationWarnings
1159 || writeExtrapolationWarnings
1160 || stopOnExtrapolationWarnings;
1161
1162 log << "*****************************************"
1163 "**************************************\n";
1164 return;
1165}
1166
1168{
1169 double rCutScreen = screeningFunction.getOuter();
1170 cutoffs.resize(numElements);
1171 for(auto const& e : elements)
1172 {
1173 e.getCutoffRadii(cutoffs.at(e.getIndex()));
1174 if ( rCutScreen > 0 &&
1175 !vectorContains(cutoffs.at(e.getIndex()), rCutScreen) )
1176 {
1177 cutoffs.at(e.getIndex()).push_back(rCutScreen);
1178 }
1179 }
1180}
1181
1183{
1184 log << "\n";
1185 log << "*** SETUP: NEURAL NETWORKS **************"
1186 "**************************************\n";
1187 log << "\n";
1188
1189 string id;
1190
1191 // Some NNP types require extra NNs.
1193 {
1194 id = "elec";
1195 nnk.push_back(id);
1196 nns[id].id = id;
1197 nns.at(id).name = "electronegativity";
1198 nns.at(id).weightFileFormat = "weightse.%03zu.data";
1199 nns.at(id).keywordSuffix = "_electrostatic";
1200 nns.at(id).keywordSuffix2 = "_charge";
1201 }
1202 else if(nnpType == NNPType::HDNNP_Q)
1203 {
1204 id = "elec";
1205 nnk.push_back(id);
1206 nns[id].id = id;
1207 nns.at(id).name = "charge";
1208 nns.at(id).weightFileFormat = "weightse.%03zu.data";
1209 nns.at(id).keywordSuffix = "_electrostatic";
1210 nns.at(id).keywordSuffix2 = "_charge";
1211 }
1212
1213 // All NNP types contain a short range NN.
1214 id = "short";
1215 nnk.push_back(id);
1216 nns[id].id = id;
1217 nns.at(id).name = "short range";
1218 nns.at(id).weightFileFormat = "weights.%03zu.data";
1219 nns.at(id).keywordSuffix = "_short";
1220 nns.at(id).keywordSuffix2 = "_short";
1221
1222 // Loop over all NNs and set global properties.
1223 for (auto& k : nnk)
1224 {
1225 // Each elements number of hidden layers.
1226 size_t globalNumHiddenLayers = 0;
1227 // Abbreviation for current NN.
1228 NNSetup& nn = nns.at(k);
1229 // Set size of NN topology vector.
1230 nn.topology.resize(numElements);
1231 // First, check for global number of hidden layers.
1232 string keyword = "global_hidden_layers" + nn.keywordSuffix;
1233 if (settings.keywordExists(keyword))
1234 {
1235 globalNumHiddenLayers = atoi(settings[keyword].c_str());
1236 for (auto& t : nn.topology)
1237 {
1238 t.numLayers = globalNumHiddenLayers + 2;
1239 }
1240 }
1241 // Now, check for per-element number of hidden layers.
1242 keyword = "element_hidden_layers" + nn.keywordSuffix;
1243 if (settings.keywordExists(keyword))
1244 {
1245 settings::Settings::KeyRange r = settings.getValues(keyword);
1246 for (settings::Settings::KeyMap::const_iterator it = r.first;
1247 it != r.second; ++it)
1248 {
1249 vector<string> args = split(reduce(it->second.first));
1250 size_t const e = elementMap[args.at(0)];
1251 size_t const n = atoi(args.at(1).c_str());
1252 nn.topology.at(e).numLayers = n + 2;
1253 }
1254 }
1255 // Check whether user has set all NN's number of layers correctly.
1256 for (auto& t : nn.topology)
1257 {
1258 if (t.numLayers == 0)
1259 {
1260 throw runtime_error("ERROR: Number of neural network hidden "
1261 "layers unset for some elements.\n");
1262 }
1263 }
1264 // Finally, allocate NN topologies data.
1265 for (auto& t : nn.topology)
1266 {
1267 t.numNeuronsPerLayer.resize(t.numLayers, 0);
1268 t.activationFunctionsPerLayer.resize(t.numLayers,
1270 }
1271
1272 // Now read global number of neurons and activation functions.
1273 vector<string> globalNumNeuronsPerHiddenLayer;
1274 keyword = "global_nodes" + nn.keywordSuffix;
1275 if (settings.keywordExists(keyword))
1276 {
1277 globalNumNeuronsPerHiddenLayer = split(reduce(settings[keyword]));
1278 if (globalNumHiddenLayers != globalNumNeuronsPerHiddenLayer.size())
1279 {
1280 throw runtime_error(strpr("ERROR: Inconsistent global NN "
1281 "topology keyword \"%s\".\n",
1282 keyword.c_str()));
1283 }
1284 }
1285 vector<string> globalActivationFunctions;
1286 keyword = "global_activation" + nn.keywordSuffix;
1287 if (settings.keywordExists(keyword))
1288 {
1289 globalActivationFunctions = split(reduce(settings[keyword]));
1290 if (globalNumHiddenLayers != globalActivationFunctions.size() - 1)
1291 {
1292 throw runtime_error(strpr("ERROR: Inconsistent global NN "
1293 "topology keyword \"%s\".\n",
1294 keyword.c_str()));
1295 }
1296 }
1297 // Set global number of neurons and activation functions if provided.
1298 bool globalNumNeurons = (globalNumNeuronsPerHiddenLayer.size() != 0);
1299 bool globalActivation = (globalActivationFunctions.size() != 0);
1300 for (size_t i = 0; i < numElements; ++i)
1301 {
1302 NNSetup::Topology& t = nn.topology.at(i);
1303 size_t const nsf = elements.at(i).numSymmetryFunctions();
1304 // Set input layer. Number of input layer neurons depends on NNP
1305 // type and NN purpose.
1307 {
1308 // Can assume NN id is "short".
1309 t.numNeuronsPerLayer.at(0) = nsf;
1310 }
1311 else if (nnpType == NNPType::HDNNP_4G ||
1313 {
1314 // NN with id "elec" requires only SFs.
1315 if (k == "elec") t.numNeuronsPerLayer.at(0) = nsf;
1316 // "short" NN needs extra charge neuron.
1317 else if (k == "short") t.numNeuronsPerLayer.at(0) = nsf + 1;
1318 }
1319 // Set dummy input neuron activation function.
1321 // Set output layer. Assume single output neuron.
1322 t.numNeuronsPerLayer.at(t.numLayers - 1) = 1;
1323 // If this element's NN does not use the global number of hidden
1324 // layers it makes no sense to set the global number of hidden
1325 // neurons or activation functions. Hence, skip the settings here,
1326 // appropriate settings should follow later in the per-element
1327 // section.
1328 if ((size_t)t.numLayers != globalNumHiddenLayers + 2) continue;
1329 for (int j = 1; j < t.numLayers; ++j)
1330 {
1331 if ((j == t.numLayers - 1) && globalActivation)
1332 {
1334 globalActivationFunctions.at(t.numLayers - 2));
1335 }
1336 else
1337 {
1338 if (globalNumNeurons)
1339 {
1340 t.numNeuronsPerLayer.at(j) = atoi(
1341 globalNumNeuronsPerHiddenLayer.at(j - 1).c_str());
1342 }
1343 if (globalActivation)
1344 {
1347 globalActivationFunctions.at(j - 1));
1348 }
1349 }
1350 }
1351 }
1352 // Override global number of neurons with per-element keyword.
1353 keyword = "element_nodes" + nn.keywordSuffix;
1354 if (settings.keywordExists(keyword))
1355 {
1356 settings::Settings::KeyRange r = settings.getValues(keyword);
1357 for (settings::Settings::KeyMap::const_iterator it = r.first;
1358 it != r.second; ++it)
1359 {
1360 vector<string> args = split(reduce(it->second.first));
1361 size_t e = elementMap[args.at(0)];
1362 size_t n = args.size() - 1;
1363 NNSetup::Topology& t = nn.topology.at(e);
1364 if ((size_t)t.numLayers != n + 2)
1365 {
1366 throw runtime_error(strpr("ERROR: Inconsistent per-element"
1367 " NN topology keyword \"%s\".\n",
1368 keyword.c_str()));
1369 }
1370 for (int j = 1; j < t.numLayers - 2; ++j)
1371 {
1372 t.numNeuronsPerLayer.at(j) = atoi(args.at(j).c_str());
1373 }
1374 }
1375 }
1376 // Override global activation functions with per-element keyword.
1377 keyword = "element_activation" + nn.keywordSuffix;
1378 if (settings.keywordExists(keyword))
1379 {
1380 settings::Settings::KeyRange r = settings.getValues(keyword);
1381 for (settings::Settings::KeyMap::const_iterator it = r.first;
1382 it != r.second; ++it)
1383 {
1384 vector<string> args = split(reduce(it->second.first));
1385 size_t e = elementMap[args.at(0)];
1386 size_t n = args.size() - 1;
1387 NNSetup::Topology& t = nn.topology.at(e);
1388 if ((size_t)t.numLayers != n + 1)
1389 {
1390 throw runtime_error(strpr("ERROR: Inconsistent per-element"
1391 " NN topology keyword \"%s\".\n",
1392 keyword.c_str()));
1393 }
1394 for (int j = 1; j < t.numLayers - 1; ++j)
1395 {
1397 activationFromString(args.at(j).c_str());
1398 }
1399 }
1400 }
1401
1402 // Finally check everything for any unset NN property.
1403 for (size_t i = 0; i < numElements; ++i)
1404 {
1405 NNSetup::Topology const& t = nn.topology.at(i);
1406 for (int j = 0; j < t.numLayers; ++j)
1407 {
1408 if (t.numNeuronsPerLayer.at(j) == 0)
1409 {
1410 throw runtime_error(strpr(
1411 "ERROR: NN \"%s\", element %2s: number of "
1412 "neurons for layer %d unset.\n",
1413 nn.id.c_str(),
1414 elements.at(i).getSymbol().c_str(),
1415 j));
1416 }
1417 if (t.activationFunctionsPerLayer.at(j)
1419 {
1420 throw runtime_error(strpr(
1421 "ERROR: NN \"%s\", element %2s: activation "
1422 "functions for layer %d unset.\n",
1423 nn.id.c_str(),
1424 elements.at(i).getSymbol().c_str(),
1425 j));
1426 }
1427 }
1428 }
1429 }
1430
1431
1432 bool normalizeNeurons = settings.keywordExists("normalize_nodes");
1433 log << strpr("Normalize neurons (all elements): %d\n",
1434 (int)normalizeNeurons);
1435 log << "-----------------------------------------"
1436 "--------------------------------------\n";
1437
1438 // Finally, allocate all neural networks.
1439 for (auto& k : nnk)
1440 {
1441 for (size_t i = 0; i < numElements; ++i)
1442 {
1443 Element& e = elements.at(i);
1444 NNSetup::Topology const& t = nns.at(k).topology.at(i);
1445 e.neuralNetworks.emplace(
1446 piecewise_construct,
1447 forward_as_tuple(k),
1448 forward_as_tuple(
1449 t.numLayers,
1450 t.numNeuronsPerLayer.data(),
1451 t.activationFunctionsPerLayer.data()));
1452 e.neuralNetworks.at(k).setNormalizeNeurons(normalizeNeurons);
1453 log << strpr("Atomic %s NN for "
1454 "element %2s :\n",
1455 nns.at(k).name.c_str(),
1456 e.getSymbol().c_str());
1457 log << e.neuralNetworks.at(k).info();
1458 log << "-----------------------------------------"
1459 "--------------------------------------\n";
1460 }
1461 }
1462
1463 log << "*****************************************"
1464 "**************************************\n";
1465
1466 return;
1467}
1468
1469void Mode::setupNeuralNetworkWeights(map<string, string> fileNameFormats)
1470{
1471 setupNeuralNetworkWeights("", fileNameFormats);
1472 return;
1473}
1474
1475void Mode::setupNeuralNetworkWeights(string directoryPrefix,
1476 map<string, string> fileNameFormats)
1477{
1478 log << "\n";
1479 log << "*** SETUP: NEURAL NETWORK WEIGHTS *******"
1480 "**************************************\n";
1481 log << "\n";
1482
1483 for (auto k : nnk)
1484 {
1485 string actualFileNameFormat;
1486 if (fileNameFormats.find(k) != fileNameFormats.end())
1487 {
1488 actualFileNameFormat = fileNameFormats.at(k);
1489 }
1490 else actualFileNameFormat = nns.at(k).weightFileFormat;
1491 actualFileNameFormat = directoryPrefix + actualFileNameFormat;
1492 log << strpr("%s weight file name format: %s\n",
1493 cap(nns.at(k).name).c_str(),
1494 actualFileNameFormat.c_str());
1495 readNeuralNetworkWeights(k, actualFileNameFormat);
1496 }
1497
1498 log << "*****************************************"
1499 "**************************************\n";
1500
1501 return;
1502}
1503
1505 bool const derivatives)
1506{
1507 // Skip calculation for whole structure if results are already saved.
1508 if (structure.hasSymmetryFunctionDerivatives) return;
1509 if (structure.hasSymmetryFunctions && !derivatives) return;
1510
1511 Atom* a = NULL;
1512 Element* e = NULL;
1513#ifdef _OPENMP
1514 #pragma omp parallel for private (a, e)
1515#endif
1516 for (size_t i = 0; i < structure.atoms.size(); ++i)
1517 {
1518 // Pointer to atom.
1519 a = &(structure.atoms.at(i));
1520
1521 // Skip calculation for individual atom if results are already saved.
1522 if (a->hasSymmetryFunctionDerivatives) continue;
1523 if (a->hasSymmetryFunctions && !derivatives) continue;
1524
1525 // Inform atom if extra charge neuron is present in short-range NN.
1526 if (nnpType == NNPType::HDNNP_4G ||
1528
1529 // Get element of atom and set number of symmetry functions.
1530 e = &(elements.at(a->element));
1532 if (derivatives)
1533 {
1536 }
1537#ifndef N2P2_NO_SF_CACHE
1539#endif
1540
1541#ifndef N2P2_NO_NEIGH_CHECK
1542 // Check if atom has low number of neighbors.
1543 size_t numNeighbors = a->calculateNumNeighbors(
1544 minCutoffRadius.at(e->getIndex()));
1545 if (numNeighbors < minNeighbors.at(e->getIndex()))
1546 {
1547 log << strpr("WARNING: Structure %6zu Atom %6zu : %zu "
1548 "neighbors.\n",
1549 a->indexStructure,
1550 a->index,
1551 numNeighbors);
1552 }
1553#endif
1554
1555 // Allocate symmetry function data vectors in atom.
1556 a->allocate(derivatives, maxCutoffRadius);
1557
1558 // Calculate symmetry functions (and derivatives).
1559 e->calculateSymmetryFunctions(*a, derivatives);
1560
1561 // Remember that symmetry functions of this atom have been calculated.
1562 a->hasSymmetryFunctions = true;
1563 if (derivatives) a->hasSymmetryFunctionDerivatives = true;
1564 }
1565
1566 // If requested, check extrapolation warnings or update statistics.
1567 // Needed to shift this out of the loop above to make it thread-safe.
1569 {
1570 for (size_t i = 0; i < structure.atoms.size(); ++i)
1571 {
1572 a = &(structure.atoms.at(i));
1573 e = &(elements.at(a->element));
1575 }
1576 }
1577
1578 // Remember that symmetry functions of this structure have been calculated.
1579 structure.hasSymmetryFunctions = true;
1580 if (derivatives) structure.hasSymmetryFunctionDerivatives = true;
1581
1582 return;
1583}
1584
1586 bool const derivatives)
1587{
1588 // Skip calculation for whole structure if results are already saved.
1589 if (structure.hasSymmetryFunctionDerivatives) return;
1590 if (structure.hasSymmetryFunctions && !derivatives) return;
1591
1592 Atom* a = NULL;
1593 Element* e = NULL;
1594#ifdef _OPENMP
1595 #pragma omp parallel for private (a, e)
1596#endif
1597 for (size_t i = 0; i < structure.atoms.size(); ++i)
1598 {
1599 // Pointer to atom.
1600 a = &(structure.atoms.at(i));
1601
1602 // Skip calculation for individual atom if results are already saved.
1603 if (a->hasSymmetryFunctionDerivatives) continue;
1604 if (a->hasSymmetryFunctions && !derivatives) continue;
1605
1606 // Inform atom if extra charge neuron is present in short-range NN.
1607 if (nnpType == NNPType::HDNNP_4G ||
1609
1610 // Get element of atom and set number of symmetry functions.
1611 e = &(elements.at(a->element));
1613 if (derivatives)
1614 {
1617 }
1618#ifndef N2P2_NO_SF_CACHE
1620#endif
1621
1622#ifndef N2P2_NO_NEIGH_CHECK
1623 // Check if atom has low number of neighbors.
1624 size_t numNeighbors = a->calculateNumNeighbors(
1625 minCutoffRadius.at(e->getIndex()));
1626 if (numNeighbors < minNeighbors.at(e->getIndex()))
1627 {
1628 log << strpr("WARNING: Structure %6zu Atom %6zu : %zu "
1629 "neighbors.\n",
1630 a->indexStructure,
1631 a->index,
1632 numNeighbors);
1633 }
1634#endif
1635
1636 // Allocate symmetry function data vectors in atom.
1637 a->allocate(derivatives, maxCutoffRadius);
1638
1639 // Calculate symmetry functions (and derivatives).
1640 e->calculateSymmetryFunctionGroups(*a, derivatives);
1641
1642 // Remember that symmetry functions of this atom have been calculated.
1643 a->hasSymmetryFunctions = true;
1644 if (derivatives) a->hasSymmetryFunctionDerivatives = true;
1645 }
1646
1647 // If requested, check extrapolation warnings or update statistics.
1648 // Needed to shift this out of the loop above to make it thread-safe.
1650 {
1651 for (size_t i = 0; i < structure.atoms.size(); ++i)
1652 {
1653 a = &(structure.atoms.at(i));
1654 e = &(elements.at(a->element));
1656 }
1657 }
1658
1659 // Remember that symmetry functions of this structure have been calculated.
1660 structure.hasSymmetryFunctions = true;
1661 if (derivatives) structure.hasSymmetryFunctionDerivatives = true;
1662
1663 return;
1664}
1665
1667 bool const derivatives,
1668 string id)
1669{
1670 if (id == "") id = nnk.front();
1671
1673 {
1674 for (vector<Atom>::iterator it = structure.atoms.begin();
1675 it != structure.atoms.end(); ++it)
1676 {
1677 NeuralNetwork& nn = elements.at(it->element)
1678 .neuralNetworks.at(id);
1679 nn.setInput(&((it->G).front()));
1680 nn.propagate();
1681 if (derivatives)
1682 {
1683 nn.calculateDEdG(&((it->dEdG).front()));
1684 }
1685 nn.getOutput(&(it->energy));
1686 }
1687 }
1688 else if (nnpType == NNPType::HDNNP_4G)
1689 {
1690 if (id == "elec")
1691 {
1692 for (auto& a : structure.atoms)
1693 {
1694 NeuralNetwork& nn = elements.at(a.element)
1695 .neuralNetworks.at(id);
1696 nn.setInput(&((a.G).front()));
1697 nn.propagate();
1698 if (derivatives)
1699 {
1700 // Calculation of dEdG is identical to dChidG
1701 nn.calculateDEdG(&((a.dChidG).front()));
1702 if (normalize)
1703 {
1704 for (auto& dChidGi : a.dChidG)
1705 dChidGi = normalized("negativity", dChidGi);
1706 }
1707 }
1708 nn.getOutput(&(a.chi));
1709 //log << strpr("Atom %5zu (%2s) chi: %24.16E\n",
1710 // a.index, elementMap[a.element].c_str(), a.chi);
1711 if (normalize) a.chi = normalized("negativity", a.chi);
1712 }
1713 }
1714 else if (id == "short")
1715 {
1716 for (auto& a : structure.atoms)
1717 {
1718 NeuralNetwork& nn = elements.at(a.element)
1719 .neuralNetworks.at(id);
1720 // TODO: This part should simplify with improved NN class.
1721 for (size_t i = 0; i < a.G.size(); ++i)
1722 {
1723 nn.setInput(i, a.G.at(i));
1724 }
1725 // Set additional charge neuron.
1726 nn.setInput(a.G.size(), a.charge);
1727 nn.propagate();
1728 if (derivatives)
1729 {
1730 // last element of vector dEdG is dEdQ
1731 nn.calculateDEdG(&((a.dEdG).front()));
1732 }
1733 nn.getOutput(&(a.energy));
1734 //log << strpr("Atom %5zu (%2s) energy: %24.16E\n",
1735 // a.index, elementMap[a.element].c_str(), a.energy);
1736 }
1737 }
1738 }
1739 else if (nnpType == NNPType::HDNNP_Q)
1740 {
1741 // Ignore ID, both NNs are computed here.
1742 for (vector<Atom>::iterator it = structure.atoms.begin();
1743 it != structure.atoms.end(); ++it)
1744 {
1745 // First the charge NN.
1746 NeuralNetwork& nnCharge = elements.at(it->element)
1747 .neuralNetworks.at("elec");
1748 nnCharge.setInput(&((it->G).front()));
1749 nnCharge.propagate();
1750 if (derivatives)
1751 {
1752 nnCharge.calculateDEdG(&((it->dQdG).front()));
1753 }
1754 nnCharge.getOutput(&(it->charge));
1755
1756 // Now the short-range NN (have to set input neurons individually).
1757 NeuralNetwork& nnShort = elements.at(it->element)
1758 .neuralNetworks.at("short");
1759 // TODO: This part should simplify with improved NN class.
1760 for (size_t i = 0; i < it->G.size(); ++i)
1761 {
1762 nnShort.setInput(i, it->G.at(i));
1763 }
1764 // Set additional charge neuron.
1765 nnShort.setInput(it->G.size(), it->charge);
1766 nnShort.propagate();
1767 if (derivatives)
1768 {
1769 nnShort.calculateDEdG(&((it->dEdG).front()));
1770 }
1771 nnShort.getOutput(&(it->energy));
1772 }
1773 }
1774
1775 return;
1776}
1777
1778// TODO: Make this const?
1780 bool const derivativesElec)
1781{
1782 Structure& s = structure;
1783
1784 // Prepare hardness vector and precalculate gamma(i, j).
1785 VectorXd hardness(numElements);
1786 MatrixXd gammaSqrt2(numElements, numElements);
1787 VectorXd sigmaSqrtPi(numElements);
1788 // In case of natural units and no normalization
1789 double fourPiEps = 1;
1790 // In case of natural units but with normalization. Other units currently
1791 // not supported.
1792 if (normalize) fourPiEps = pow(convCharge, 2) / (convLength * convEnergy);
1793
1794 fourPiEps = 1;
1795
1796 for (size_t i = 0; i < numElements; ++i)
1797 {
1798 hardness(i) = elements.at(i).getHardness();
1799 double const iSigma = elements.at(i).getQsigma();
1800 sigmaSqrtPi(i) = sqrt(M_PI) * iSigma;
1801 for (size_t j = 0; j < numElements; ++j)
1802 {
1803 double const jSigma = elements.at(j).getQsigma();
1804 gammaSqrt2(i, j) = sqrt(2.0 * (iSigma * iSigma + jSigma * jSigma));
1805 }
1806 }
1807
1809 hardness,
1810 gammaSqrt2,
1811 sigmaSqrtPi,
1813 fourPiEps,
1814 erfcBuf);
1815
1816 // TODO: leave these 2 functions here or move it to e.g. forces? Needs to be
1817 // executed after calculateElectrostaticEnergy.
1818 if (derivativesElec)
1819 {
1820 s.calculateDAdrQ(ewaldSetup, gammaSqrt2, fourPiEps, erfcBuf);
1822 gammaSqrt2,
1823 sigmaSqrtPi,
1825 fourPiEps);
1826 }
1827
1828 return;
1829}
1830
1831void Mode::calculateEnergy(Structure& structure) const
1832{
1833 // Loop over all atoms and add atomic contributions to total energy.
1834 structure.energy = 0.0;
1835 structure.energyShort = 0.0;
1836 for (vector<Atom>::iterator it = structure.atoms.begin();
1837 it != structure.atoms.end(); ++it)
1838 {
1839 structure.energyShort += it->energy;
1840 }
1841 structure.energy = structure.energyShort + structure.energyElec;
1842
1843 //cout << strpr("Electrostatic energy: %24.16E\n", structure.energyElec);
1844 //cout << strpr("Short-range energy: %24.16E\n", structure.energyShort);
1845 //cout << strpr("Sum energy: %24.16E\n", structure.energy);
1846 //cout << strpr("Offset energy: %24.16E\n", getEnergyOffset(structure));
1847 //cout << "---------------------\n";
1848 //cout << strpr("Total energy: %24.16E\n", structure.energy + getEnergyOffset(structure));
1849 //cout << strpr("Reference energy: %24.16E\n", structure.energyRef + getEnergyOffset(structure));
1850 //cout << "---------------------\n";
1851 //cout << "without offset: \n";
1852 //cout << strpr("Total energy: %24.16E\n", structure.energy);
1853 //cout << strpr("Reference energy: %24.16E\n", structure.energyRef);
1854
1855 return;
1856}
1857
1858void Mode::calculateCharge(Structure& structure) const
1859{
1860 // Loop over all atoms and add atomic charge contributions to total charge.
1861 structure.charge = 0.0;
1862 for (vector<Atom>::iterator it = structure.atoms.begin();
1863 it != structure.atoms.end(); ++it)
1864 {
1865 //log << strpr("Atom %5zu (%2s) q: %24.16E\n",
1866 // it->index, elementMap[it->element].c_str(), it->charge);
1867 structure.charge += it->charge;
1868 }
1869
1870 //log << "---------------------\n";
1871 //log << strpr("Total charge: %24.16E\n", structure.charge);
1872 //log << strpr("Reference charge: %24.16E\n", structure.chargeRef);
1873
1874 return;
1875}
1876
1877void Mode::calculateForces(Structure& structure) const
1878{
1880 {
1881 throw runtime_error("WARNING: Forces are not implemented yet.\n");
1882 return;
1883 }
1884
1885 // Loop over all atoms, center atom i (ai).
1886#ifdef _OPENMP
1887 #pragma omp parallel
1888 {
1889 #pragma omp for
1890#endif
1891 for (size_t i = 0; i < structure.atoms.size(); ++i) {
1892 // Set pointer to atom.
1893 Atom &ai = structure.atoms.at(i);
1894
1895 // Reset forces.
1896 ai.f = Vec3D{};
1897
1898 // First add force contributions from atom i itself (gradient of
1899 // atomic energy E_i).
1900 ai.f += ai.calculateSelfForceShort();
1901
1902 // Now loop over all neighbor atoms j of atom i. These may hold
1903 // non-zero derivatives of their symmetry functions with respect to
1904 // atom i's coordinates. Some atoms may appear multiple times in the
1905 // neighbor list because of periodic boundary conditions. To avoid
1906 // that the same contributions are added multiple times use the
1907 // "unique neighbor" list. This list contains also the central atom
1908 // index as first entry and hence also adds contributions of periodic
1909 // images of the central atom (happens when cutoff radii larger than
1910 // cell vector lengths are used, but this is already considered in the
1911 // self-interaction).
1912 for (vector<size_t>::const_iterator it = ai.neighborsUnique.begin() + 1;
1913 it != ai.neighborsUnique.end(); ++it)
1914 {
1915 // Define shortcut for atom j (aj).
1916 Atom &aj = structure.atoms.at(*it);
1917#ifndef N2P2_FULL_SFD_MEMORY
1918 vector<vector<size_t>> const& tableFull
1919 = elements.at(aj.element).getSymmetryFunctionTable();
1920#endif
1921 // Loop over atom j's neighbors (n), atom i should be one of them.
1922 // TODO: Could implement maxCutoffRadiusSymFunc for each element and
1923 // use this instead of maxCutoffRadius.
1924 size_t const numNeighbors =
1926 for (size_t k = 0; k < numNeighbors; ++k) {
1927 Atom::Neighbor const &n = aj.neighbors[k];
1928 // If atom j's neighbor is atom i add force contributions.
1929 if (n.index == ai.index)
1930 {
1931#ifndef N2P2_FULL_SFD_MEMORY
1932 ai.f += aj.calculatePairForceShort(n, &tableFull);
1933#else
1934 ai.f += aj.calculatePairForceShort(n);
1935#endif
1936 }
1937 }
1938 }
1939 }
1940
1942 {
1943 Structure &s = structure;
1944
1945 VectorXd lambdaTotal = s.calculateForceLambdaTotal();
1946 VectorXd lambdaElec = s.calculateForceLambdaElec();
1947
1948#ifdef _OPENMP
1949 #pragma omp for
1950#endif
1951 // OpenMP 4.0 doesn't support range based loops
1952 for (size_t i = 0; i < s.numAtoms; ++i)
1953 {
1954 auto &ai = s.atoms[i];
1955 ai.f -= ai.pEelecpr;
1956 ai.fElec = -ai.pEelecpr;
1957
1958 for (size_t j = 0; j < s.numAtoms; ++j)
1959 {
1960 Atom const &aj = s.atoms.at(j);
1961
1962#ifndef N2P2_FULL_SFD_MEMORY
1963 vector<vector<size_t> > const &tableFull
1964 = elements.at(aj.element).getSymmetryFunctionTable();
1965 Vec3D dChidr = aj.calculateDChidr(ai.index,
1967 &tableFull);
1968#else
1969 Vec3D dChidr = aj.calculateDChidr(ai.index,
1971#endif
1972 ai.f -= lambdaTotal(j) * (ai.dAdrQ[j] + dChidr);
1973 ai.fElec -= lambdaElec(j) * (ai.dAdrQ[j] + dChidr);
1974 }
1975 }
1976 }
1977#ifdef _OPENMP
1978 }
1979#endif
1980 return;
1981}
1982
1983
1984void Mode::evaluateNNP(Structure& structure, bool useForces, bool useDEdG)
1985{
1986 useDEdG = (useForces || useDEdG);
1988 {
1990 ewaldSetup,
1991 screeningFunction.getOuter(),
1994 }
1995 // TODO: For the moment sort neighbors only for 4G-HDNNPs (breaks some
1996 // CI tests because of small numeric changes).
1997 else structure.calculateNeighborList(maxCutoffRadius, false);
1998
1999#ifdef N2P2_NO_SF_GROUPS
2000 calculateSymmetryFunctions(structure, true);
2001#else
2002 calculateSymmetryFunctionGroups(structure, useForces);
2003#endif
2004
2005 // Needed if useForces == false but sensitivity data is requested.
2006 if (useDEdG && !useForces)
2007 {
2008 // Manually allocate dEdG vectors.
2009 for (auto& a : structure.atoms)
2010 {
2012 {
2013 a.dEdG.resize(a.numSymmetryFunctions + 1, 0.0);
2014 a.dChidG.resize(a.numSymmetryFunctions, 0.0);
2015 }
2016 else
2017 a.dEdG.resize(a.numSymmetryFunctions, 0.0);
2018 }
2019 }
2020
2021 calculateAtomicNeuralNetworks(structure, useDEdG);
2023 {
2024 chargeEquilibration(structure, useForces);
2025 calculateAtomicNeuralNetworks(structure, useDEdG, "short");
2026 }
2027 calculateEnergy(structure);
2028 if (nnpType == NNPType::HDNNP_4G ||
2030 calculateCharge(structure);
2031 if (useForces) calculateForces(structure);
2032}
2033
2034void Mode::addEnergyOffset(Structure& structure, bool ref)
2035{
2036 for (size_t i = 0; i < numElements; ++i)
2037 {
2038 if (ref)
2039 {
2040 structure.energyRef += structure.numAtomsPerElement.at(i)
2041 * elements.at(i).getAtomicEnergyOffset();
2042 }
2043 else
2044 {
2045 structure.energy += structure.numAtomsPerElement.at(i)
2046 * elements.at(i).getAtomicEnergyOffset();
2047 }
2048 }
2049
2050 return;
2051}
2052
2053void Mode::removeEnergyOffset(Structure& structure, bool ref)
2054{
2055 for (size_t i = 0; i < numElements; ++i)
2056 {
2057 if (ref)
2058 {
2059 structure.energyRef -= structure.numAtomsPerElement.at(i)
2060 * elements.at(i).getAtomicEnergyOffset();
2061 }
2062 else
2063 {
2064 structure.energy -= structure.numAtomsPerElement.at(i)
2065 * elements.at(i).getAtomicEnergyOffset();
2066 }
2067 }
2068
2069 return;
2070}
2071
2072double Mode::getEnergyOffset(Structure const& structure) const
2073{
2074 double result = 0.0;
2075
2076 for (size_t i = 0; i < numElements; ++i)
2077 {
2078 result += structure.numAtomsPerElement.at(i)
2079 * elements.at(i).getAtomicEnergyOffset();
2080 }
2081
2082 return result;
2083}
2084
2085double Mode::getEnergyWithOffset(Structure const& structure, bool ref) const
2086{
2087 double result;
2088 if (ref) result = structure.energyRef;
2089 else result = structure.energy;
2090
2091 for (size_t i = 0; i < numElements; ++i)
2092 {
2093 result += structure.numAtomsPerElement.at(i)
2094 * elements.at(i).getAtomicEnergyOffset();
2095 }
2096
2097 return result;
2098}
2099
2100double Mode::normalized(string const& property, double value) const
2101{
2102 if (property == "energy") return value * convEnergy;
2103 else if (property == "force") return value * convEnergy / convLength;
2104 else if (property == "charge") return value * convCharge;
2105 else if (property == "hardness")
2106 return value * convEnergy / pow(convCharge, 2);
2107 else if (property == "negativity") return value * convEnergy / convCharge;
2108 else throw runtime_error("ERROR: Unknown property to convert to "
2109 "normalized units.\n");
2110}
2111
2112double Mode::normalizedEnergy(Structure const& structure, bool ref) const
2113{
2114 if (ref)
2115 {
2116 return (structure.energyRef - structure.numAtoms * meanEnergy)
2117 * convEnergy;
2118 }
2119 else
2120 {
2121 return (structure.energy - structure.numAtoms * meanEnergy)
2122 * convEnergy;
2123 }
2124}
2125
2126double Mode::physical(string const& property, double value) const
2127{
2128 if (property == "energy") return value / convEnergy;
2129 else if (property == "force") return value * convLength / convEnergy;
2130 else if (property == "charge") return value / convCharge;
2131 else if (property == "hardness")
2132 return value / (convEnergy / pow(convCharge, 2));
2133 else if (property == "negativity") return value * convCharge / convEnergy;
2134 else throw runtime_error("ERROR: Unknown property to convert to physical "
2135 "units.\n");
2136}
2137
2138double Mode::physicalEnergy(Structure const& structure, bool ref) const
2139{
2140 if (ref)
2141 {
2142 return structure.energyRef / convEnergy + structure.numAtoms
2143 * meanEnergy;
2144 }
2145 else
2146 {
2147 return structure.energy / convEnergy + structure.numAtoms * meanEnergy;
2148 }
2149}
2150
2152{
2154
2155 return;
2156}
2157
2159{
2161
2162 return;
2163}
2164
2166{
2167 if (nnpType != NNPType::HDNNP_4G) return;
2168 ewaldSetup.logEwaldCutoffs(log, convLength);
2169}
2170
2172{
2173 for (vector<Element>::iterator it = elements.begin();
2174 it != elements.end(); ++it)
2175 {
2176 it->statistics.resetExtrapolationWarnings();
2177 }
2178
2179 return;
2180}
2181
2183{
2184 size_t numExtrapolationWarnings = 0;
2185
2186 for (vector<Element>::const_iterator it = elements.begin();
2187 it != elements.end(); ++it)
2188 {
2189 numExtrapolationWarnings +=
2190 it->statistics.countExtrapolationWarnings();
2191 }
2192
2193 return numExtrapolationWarnings;
2194}
2195
2196vector<size_t> Mode::getNumSymmetryFunctions() const
2197{
2198 vector<size_t> v;
2199
2200 for (vector<Element>::const_iterator it = elements.begin();
2201 it != elements.end(); ++it)
2202 {
2203 v.push_back(it->numSymmetryFunctions());
2204 }
2205
2206 return v;
2207}
2208
2209bool Mode::settingsKeywordExists(std::string const& keyword) const
2210{
2211 return settings.keywordExists(keyword);
2212}
2213
2214string Mode::settingsGetValue(std::string const& keyword) const
2215{
2216 return settings.getValue(keyword);
2217}
2218
2219
2220void Mode::writePrunedSettingsFile(vector<size_t> prune, string fileName) const
2221{
2222 ofstream file(fileName.c_str());
2223 vector<string> settingsLines = settings.getSettingsLines();
2224 for (size_t i = 0; i < settingsLines.size(); ++i)
2225 {
2226 if (find(prune.begin(), prune.end(), i) != prune.end())
2227 {
2228 file << "# ";
2229 }
2230 file << settingsLines.at(i) << '\n';
2231 }
2232 file.close();
2233
2234 return;
2235}
2236
2237void Mode::writeSettingsFile(ofstream* const& file) const
2238{
2239 settings.writeSettingsFile(file);
2240
2241 return;
2242}
2243
2244vector<size_t> Mode::pruneSymmetryFunctionsRange(double threshold)
2245{
2246 vector<size_t> prune;
2247
2248 // Check if symmetry functions have low range.
2249 for (vector<Element>::const_iterator it = elements.begin();
2250 it != elements.end(); ++it)
2251 {
2252 for (size_t i = 0; i < it->numSymmetryFunctions(); ++i)
2253 {
2254 SymFnc const& s = it->getSymmetryFunction(i);
2255 if (fabs(s.getGmax() - s.getGmin()) < threshold)
2256 {
2257 prune.push_back(it->getSymmetryFunction(i).getLineNumber());
2258 }
2259 }
2260 }
2261
2262 return prune;
2263}
2264
2266 double threshold,
2267 vector<vector<double> > sensitivity)
2268{
2269 vector<size_t> prune;
2270
2271 for (size_t i = 0; i < numElements; ++i)
2272 {
2273 for (size_t j = 0; j < elements.at(i).numSymmetryFunctions(); ++j)
2274 {
2275 if (sensitivity.at(i).at(j) < threshold)
2276 {
2277 prune.push_back(
2278 elements.at(i).getSymmetryFunction(j).getLineNumber());
2279 }
2280 }
2281 }
2282
2283 return prune;
2284}
2285
2287 string const& fileNameFormat)
2288{
2289 for (vector<Element>::iterator it = elements.begin();
2290 it != elements.end(); ++it)
2291 {
2292 string fileName = strpr(fileNameFormat.c_str(),
2293 it->getAtomicNumber());
2294 log << strpr("Setting weights for element %2s from file: %s\n",
2295 it->getSymbol().c_str(),
2296 fileName.c_str());
2297 vector<double> weights = readColumnsFromFile(fileName,
2298 vector<size_t>(1, 0)
2299 ).at(0);
2300 NeuralNetwork& nn = it->neuralNetworks.at(id);
2301 nn.setConnections(&(weights.front()));
2302 }
2303
2304 return;
2305}
CutoffType
List of available cutoff function types.
Contains element-specific data.
Definition Element.h:39
void calculateSymmetryFunctions(Atom &atom, bool const derivatives) const
Calculate symmetry functions.
Definition Element.cpp:436
void calculateSymmetryFunctionGroups(Atom &atom, bool const derivatives) const
Calculate symmetry functions via groups.
Definition Element.cpp:449
std::vector< std::size_t > getCacheSizes() const
Get cache sizes for each neighbor atom element.
Definition Element.cpp:563
std::size_t getIndex() const
Get index.
Definition Element.h:305
void setCacheIndices(std::vector< std::vector< SFCacheList > > cacheLists)
Set cache indices for all symmetry functions of this element.
Definition Element.cpp:544
std::size_t updateSymmetryFunctionStatistics(Atom const &atom)
Update symmetry function statistics.
Definition Element.cpp:462
SymFnc const & getSymmetryFunction(std::size_t index) const
Get symmetry function instance.
Definition Element.h:358
std::string getSymbol() const
Get symbol.
Definition Element.h:330
std::map< std::string, NeuralNetwork > neuralNetworks
Neural networks for this element.
Definition Element.h:247
std::size_t numSymmetryFunctions() const
Get number of symmetry functions.
Definition Element.h:353
std::vector< std::size_t > const & getSymmetryFunctionNumTable() const
Get number of relevant symmetry functions per element.
Definition Element.h:336
Log log
Global log file.
Definition Mode.h:629
double physicalEnergy(Structure const &structure, bool ref=true) const
Undo normalization for a given energy of structure.
Definition Mode.cpp:2138
bool checkExtrapolationWarnings
Definition Mode.h:666
std::vector< std::vector< double > > cutoffs
Matrix storing all symmetry function cut-offs for all elements.
Definition Mode.h:689
bool normalize
Definition Mode.h:665
double convEnergy
Definition Mode.h:673
void setupNormalization(bool standalone=true)
Set up normalization.
Definition Mode.cpp:239
ElementMap elementMap
Global element map, populated by setupElementMap().
Definition Mode.h:627
NNPType nnpType
Definition Mode.h:664
double fourPiEps
Definition Mode.h:676
std::vector< double > minCutoffRadius
Definition Mode.h:669
void addEnergyOffset(Structure &structure, bool ref=true)
Add atomic energy offsets to reference energy.
Definition Mode.cpp:2034
void initialize()
Write welcome message with version information.
Definition Mode.cpp:56
double convLength
Definition Mode.h:674
void logEwaldCutoffs()
Logs Ewald params whenever they change.
Definition Mode.cpp:2165
virtual void setupElementMap()
Set up the element map.
Definition Mode.cpp:302
void readNeuralNetworkWeights(std::string const &id, std::string const &fileName)
Read in weights for a specific type of neural network.
Definition Mode.cpp:2286
virtual void setupNeuralNetwork()
Set up neural networks for all elements.
Definition Mode.cpp:1182
virtual void setupSymmetryFunctionCache(bool verbose=false)
Set up symmetry function cache.
Definition Mode.cpp:1008
double maxCutoffRadius
Definition Mode.h:670
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:213
NNPType
Definition Mode.h:91
@ HDNNP_2G
Short range NNP (2G-HDNNP).
Definition Mode.h:93
@ HDNNP_Q
Short range NNP with charge NN, no electrostatics/Qeq (M.
Definition Mode.h:103
@ HDNNP_4G
NNP with electrostatics and non-local charge transfer (4G-HDNNP).
Definition Mode.h:95
std::vector< std::string > nnk
Definition Mode.h:684
virtual void setupElectrostatics(bool initialHardness=false, std::string directoryPrefix="", std::string fileNameFormat="hardness.%03zu.data")
Set up electrostatics related stuff (hardness, screening, ...).
Definition Mode.cpp:369
std::vector< std::size_t > pruneSymmetryFunctionsRange(double threshold)
Prune symmetry functions according to their range and write settings file.
Definition Mode.cpp:2244
double meanEnergy
Definition Mode.h:672
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:1469
ScreeningFunction screeningFunction
Definition Mode.h:682
SymFnc::ScalingType scalingType
Definition Mode.h:680
void calculateAtomicNeuralNetworks(Structure &structure, bool const derivatives, std::string id="")
Calculate a single atomic neural network for a given atom and nn type.
Definition Mode.cpp:1666
double physical(std::string const &property, double value) const
Undo normalization for a given property.
Definition Mode.cpp:2126
void writePrunedSettingsFile(std::vector< std::size_t > prune, std::string fileName="output.nn") const
Copy settings file but comment out lines provided.
Definition Mode.cpp:2220
virtual void setupSymmetryFunctionGroups()
Set up symmetry function groups.
Definition Mode.cpp:866
std::vector< Element > elements
Definition Mode.h:683
double normalized(std::string const &property, double value) const
Apply normalization to given property.
Definition Mode.cpp:2100
std::size_t numElements
Definition Mode.h:667
void convertToNormalizedUnits(Structure &structure) const
Convert one structure to normalized units.
Definition Mode.cpp:2151
void calculateEnergy(Structure &structure) const
Calculate potential energy for a given structure.
Definition Mode.cpp:1831
void calculateSymmetryFunctionGroups(Structure &structure, bool const derivatives)
Calculate all symmetry function groups for all atoms in given structure.
Definition Mode.cpp:1585
settings::Settings settings
Definition Mode.h:679
void setupSymmetryFunctionScalingNone()
Set up "empy" symmetry function scaling.
Definition Mode.cpp:716
std::string settingsGetValue(std::string const &keyword) const
Get value for given keyword in Settings instance.
Definition Mode.cpp:2214
virtual void setupSymmetryFunctionScaling(std::string const &fileName="scaling.data")
Set up symmetry function scaling from file.
Definition Mode.cpp:736
double convCharge
Definition Mode.h:675
std::vector< std::size_t > getNumSymmetryFunctions() const
Get number of symmetry functions per element.
Definition Mode.cpp:2196
ErfcBuf erfcBuf
Definition Mode.h:690
std::size_t getNumExtrapolationWarnings() const
Count total number of extrapolation warnings encountered for all elements and symmetry functions.
Definition Mode.cpp:2182
void removeEnergyOffset(Structure &structure, bool ref=true)
Remove atomic energy offsets from reference energy.
Definition Mode.cpp:2053
bool settingsKeywordExists(std::string const &keyword) const
Check if keyword was found in settings file.
Definition Mode.cpp:2209
Log log
Global log file.
Definition Mode.h:629
KspaceGrid kspaceGrid
Definition Mode.h:678
void resetExtrapolationWarnings()
Erase all extrapolation warnings and reset counters.
Definition Mode.cpp:2171
double getEnergyWithOffset(Structure const &structure, bool ref=true) const
Add atomic energy offsets and return energy.
Definition Mode.cpp:2085
std::vector< std::size_t > minNeighbors
Definition Mode.h:668
void convertToPhysicalUnits(Structure &structure) const
Convert one structure to physical units.
Definition Mode.cpp:2158
void calculateCharge(Structure &structure) const
Calculate total charge for a given structure.
Definition Mode.cpp:1858
virtual void setupSymmetryFunctions()
Set up all symmetry functions.
Definition Mode.cpp:639
void chargeEquilibration(Structure &structure, bool const derivativesElec)
Perform global charge equilibration method.
Definition Mode.cpp:1779
void calculateSymmetryFunctions(Structure &structure, bool const derivatives)
Calculate all symmetry functions for all atoms in given structure.
Definition Mode.cpp:1504
void evaluateNNP(Structure &structure, bool useForces=true, bool useDEdG=true)
Evaluate neural network potential (includes total energy, optionally forces and in some cases charges...
Definition Mode.cpp:1984
double cutoffAlpha
Definition Mode.h:671
void setupCutoffMatrix()
Setup matrix storing all symmetry function cut-offs for each element.
Definition Mode.cpp:1167
std::map< std::string, NNSetup > nns
Definition Mode.h:686
EwaldSetup ewaldSetup
Definition Mode.h:677
void calculateForces(Structure &structure) const
Calculate forces for all atoms in given structure.
Definition Mode.cpp:1877
void setupSymmetryFunctionStatistics(bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
Set up symmetry function statistics collection.
Definition Mode.cpp:1127
void setupCutoff()
Set up cutoff function for all symmetry functions.
Definition Mode.cpp:544
void loadSettingsFile(std::string const &fileName="input.nn")
Open settings file and load all keywords into memory.
Definition Mode.cpp:162
std::vector< std::size_t > pruneSymmetryFunctionsSensitivity(double threshold, std::vector< std::vector< double > > sensitivity)
Prune symmetry functions with sensitivity analysis data.
Definition Mode.cpp:2265
double getEnergyOffset(Structure const &structure) const
Get atomic energy offset for given structure.
Definition Mode.cpp:2072
void setupSymmetryFunctionMemory(bool verbose=false)
Extract required memory dimensions for symmetry function derivatives.
Definition Mode.cpp:918
void writeSettingsFile(std::ofstream *const &file) const
Write complete settings file.
Definition Mode.cpp:2237
double normalizedEnergy(Structure const &structure, bool ref=true) const
Apply normalization to given energy of structure.
Definition Mode.cpp:2112
CutoffFunction::CutoffType cutoffType
Definition Mode.h:681
virtual void setupElements()
Set up all Element instances.
Definition Mode.cpp:323
This class implements a feed-forward neural network.
@ AF_UNSET
Unset activation function.
void setInput(double const *const &input) const
Set neural network input layer node values.
void setConnections(double const *const &connections)
Set neural network weights and biases.
void propagate()
Propagate input information through all layers.
void calculateDEdG(double *dEdG) const
Calculate derivative of output neuron with respect to input neurons.
void getOutput(double *output) const
Get neural network output layer node values.
Symmetry function base class.
Definition SymFnc.h:40
@ ST_SCALESIGMA
Definition SymFnc.h:64
@ ST_SCALECENTER
Definition SymFnc.h:60
@ ST_CENTER
Definition SymFnc.h:55
virtual std::string parameterLine() const =0
Give symmetry function parameters in one line.
double getGmin() const
Get private Gmin member variable.
Definition SymFnc.h:361
double getGmax() const
Get private Gmax member variable.
Definition SymFnc.h:362
std::vector< std::size_t > getIndexPerElement() const
Get private indexPerElement member variable.
Definition SymFnc.h:383
virtual std::vector< std::string > getCacheIdentifiers() const
Get unique cache identifiers.
Definition SymFnc.cpp:101
std::size_t getIndex() const
Get private index member variable.
Definition SymFnc.h:357
std::pair< KeyMap::const_iterator, KeyMap::const_iterator > KeyRange
Definition Settings.h:44
Definition Atom.h:29
@ KOLAFA_PERRAM
Method 1: Optimized in n2p2 (DOI: 10.1080/08927029208049126).
Definition EwaldSetup.h:33
@ JACKSON_CATLOW
Method 0: Used by RuNNer (DOI: 10.1080/08927028808080944).
Definition EwaldSetup.h:31
string cap(string word)
Capitalize first letter of word.
Definition utility.cpp:104
string strpr(const char *format,...)
String version of printf function.
Definition utility.cpp:90
bool vectorContains(std::vector< T > const &stdVec, T value)
Test if vector contains specified value.
Definition utility.h:166
vector< string > split(string const &input, char delimiter)
Split string at each delimiter.
Definition utility.cpp:33
map< size_t, vector< double > > readColumnsFromFile(string fileName, vector< size_t > columns, char comment)
Definition utility.cpp:247
NeuralNetwork::ActivationFunction activationFromString(std::string c)
Convert string to activation function.
@ EWALD_SUM_LAMMPS
Solver 2: Ewald summation in LAMMPS.
Definition Kspace.h:32
@ EWALD_SUM
Solver 0: Ewald summation.
Definition Kspace.h:28
@ PPPM
Solver 1: PPPM.
Definition Kspace.h:30
string reduce(string const &line, string const &whitespace, string const &fill)
Replace multiple whitespaces with fill.
Definition utility.cpp:60
Struct to store information on neighbor atoms.
Definition Atom.h:36
std::size_t index
Index of neighbor atom.
Definition Atom.h:38
Storage for a single atom.
Definition Atom.h:33
std::vector< Neighbor > neighbors
Neighbor array (maximum number defined in macros.h.
Definition Atom.h:170
std::size_t numSymmetryFunctions
Number of symmetry functions used to describe the atom environment.
Definition Atom.h:113
Vec3D calculateDChidr(size_t const atomIndexOfR, double const maxCutoffRadius, std::vector< std::vector< size_t > > const *const tableFull=nullptr) const
Calculate dChi/dr of this atom's Chi with respect to the coordinates of the given atom.
Definition Atom.cpp:403
Vec3D f
Force vector calculated by neural network.
Definition Atom.h:127
Vec3D calculatePairForceShort(Neighbor const &neighbor, std::vector< std::vector< std::size_t > > const *const tableFull=nullptr) const
Calculate force resulting from gradient of this atom's (short-ranged) energy contribution with respec...
Definition Atom.cpp:381
bool hasSymmetryFunctionDerivatives
If symmetry function derivatives are saved for this atom.
Definition Atom.h:97
std::size_t index
Index number of this atom.
Definition Atom.h:101
Vec3D calculateSelfForceShort() const
Calculate force resulting from gradient of this atom's (short-ranged) energy contribution with respec...
Definition Atom.cpp:371
std::size_t getStoredMinNumNeighbors(double const cutoffRadius) const
Return needed number of neighbors for a given cutoff radius from neighborCutoffs map.
Definition Atom.cpp:329
std::vector< std::size_t > numSymmetryFunctionDerivatives
Number of neighbor atom symmetry function derivatives per element.
Definition Atom.h:140
bool useChargeNeuron
If an additional charge neuron in the short-range NN is present.
Definition Atom.h:99
std::size_t indexStructure
Index number of structure this atom belongs to.
Definition Atom.h:103
std::size_t element
Element index of this atom.
Definition Atom.h:107
void allocate(bool all, double const maxCutoffRadius=0.0)
Allocate vectors related to symmetry functions (G, dEdG).
Definition Atom.cpp:168
bool hasSymmetryFunctions
If symmetry function values are saved for this atom.
Definition Atom.h:95
std::size_t calculateNumNeighbors(double const cutoffRadius) const
Calculate number of neighbors for a given cutoff radius.
Definition Atom.cpp:307
std::vector< std::size_t > cacheSizePerElement
Cache size for each element.
Definition Atom.h:143
std::vector< std::size_t > neighborsUnique
List of unique neighbor indices (don't count multiple PBC images).
Definition Atom.h:136
List of symmetry functions corresponding to one cache identifier.
Definition Element.h:45
std::vector< NeuralNetwork::ActivationFunction > activationFunctionsPerLayer
Activation function type per layer.
Definition Mode.h:643
std::vector< int > numNeuronsPerLayer
Number of neurons per layer.
Definition Mode.h:640
int numLayers
Number of NN layers (including input and output layer).
Definition Mode.h:638
Setup data for one neural network.
Definition Mode.h:634
std::string keywordSuffix
Suffix for keywords (NN topology related).
Definition Mode.h:656
std::vector< Topology > topology
Per-element NN topology.
Definition Mode.h:660
std::string id
NN identifier, e.g. "short", "charge",...
Definition Mode.h:650
Storage for one atomic configuration.
Definition Structure.h:39
void toPhysicalUnits(double meanEnergy, double convEnergy, double convLength, double convCharge)
Switch to physical units, shift energy and change energy, length and charge unit.
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
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.
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.
double energyShort
Short-range part of the potential energy predicted by NNP.
Definition Structure.h:88
bool hasSymmetryFunctionDerivatives
If symmetry function derivatives are saved for each atom.
Definition Structure.h:71
Eigen::VectorXd const calculateForceLambdaTotal() const
Calculate lambda_total vector which is needed for the total force calculation in 4G NN.
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
double energyRef
Reference potential energy.
Definition Structure.h:85
void toNormalizedUnits(double meanEnergy, double convEnergy, double convLength, double convCharge)
Normalize structure, shift energy and change energy, length and charge unit.
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...
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
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
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
#define N2P2_GIT_BRANCH
Definition version.h:23
#define N2P2_GIT_VERSION
Definition version.h:21
#define N2P2_VERSION
Definition version.h:20
#define N2P2_GIT_REV
Definition version.h:22