n2p2 - A neural network potential package
Loading...
Searching...
No Matches
Atom.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 "Atom.h"
18#include "utility.h"
19#include <cinttypes> // PRId64
20#include <cmath> // fabs
21#include <limits> // std::numeric_limits
22#include <stdexcept> // std::range_error
23
24using namespace std;
25using namespace nnp;
26
31 useChargeNeuron (false),
32 index (0 ),
33 indexStructure (0 ),
34 tag (0 ),
35 element (0 ),
36 numNeighbors (0 ),
39 energy (0.0 ),
40 dEelecdQ (0.0 ),
41 chi (0.0 ),
42 charge (0.0 ),
43 chargeRef (0.0 )
44{
45}
46
47#ifdef N2P2_FULL_SFD_MEMORY
48void Atom::collectDGdxia(size_t indexAtom,
49 size_t indexComponent,
50 double maxCutoffRadius)
51{
52 for (size_t i = 0; i < dGdxia.size(); i++)
53 {
54 dGdxia[i] = 0.0;
55 }
56 for (size_t i = 0; i < getStoredMinNumNeighbors(maxCutoffRadius); i++)
57 {
58 if (neighbors[i].index == indexAtom)
59 {
60 for (size_t j = 0; j < numSymmetryFunctions; ++j)
61 {
62 dGdxia[j] += neighbors[i].dGdr[j][indexComponent];
63 }
64 }
65 }
66 if (index == indexAtom)
67 {
68 for (size_t i = 0; i < numSymmetryFunctions; ++i)
69 {
70 dGdxia[i] += dGdr[i][indexComponent];
71 }
72 }
73
74 return;
75}
76#endif
77
78void Atom::toNormalizedUnits(double convEnergy,
79 double convLength,
80 double convCharge)
81{
82 energy *= convEnergy;
83 r *= convLength;
84 f *= convEnergy / convLength;
85 fRef *= convEnergy / convLength;
86 charge *= convCharge;
87 chargeRef *= convCharge;
89 {
90 for (size_t i = 0; i < numSymmetryFunctions; ++i)
91 {
92 dEdG.at(i) *= convEnergy;
93 dGdr.at(i) /= convLength;
94#ifdef N2P2_FULL_SFD_MEMORY
95 dGdxia.at(i) /= convLength;
96#endif
97 }
98 // Take care of extra charge neuron.
100 dEdG.at(numSymmetryFunctions) *= convEnergy / convCharge;
101 }
102
103 if (hasNeighborList)
104 {
105 for (vector<Neighbor>::iterator it = neighbors.begin();
106 it != neighbors.end(); ++it)
107 {
108 it->d *= convLength;
109 it->dr *= convLength;
111 {
112 for (size_t i = 0; i < dGdr.size(); ++i)
113 {
114 dGdr.at(i) /= convLength;
115 }
116 }
117 }
118 }
119
120 return;
121}
122
123void Atom::toPhysicalUnits(double convEnergy,
124 double convLength,
125 double convCharge)
126{
127 energy /= convEnergy;
128 r /= convLength;
129 f /= convEnergy / convLength;
130 fRef /= convEnergy / convLength;
131 charge /= convCharge;
132 chargeRef /= convCharge;
134 {
135 for (size_t i = 0; i < numSymmetryFunctions; ++i)
136 {
137 dEdG.at(i) /= convEnergy;
138 dGdr.at(i) *= convLength;
139#ifdef N2P2_FULL_SFD_MEMORY
140 dGdxia.at(i) *= convLength;
141#endif
142 }
143 // Take care of extra charge neuron.
144 if (useChargeNeuron)
145 dEdG.at(numSymmetryFunctions) /= convEnergy / convCharge;
146 }
147
148 if (hasNeighborList)
149 {
150 for (vector<Neighbor>::iterator it = neighbors.begin();
151 it != neighbors.end(); ++it)
152 {
153 it->d /= convLength;
154 it->dr /= convLength;
156 {
157 for (size_t i = 0; i < dGdr.size(); ++i)
158 {
159 dGdr.at(i) *= convLength;
160 }
161 }
162 }
163 }
164
165 return;
166}
167
168void Atom::allocate(bool all, double const maxCutoffRadius)
169{
170 if (numSymmetryFunctions == 0)
171 {
172 throw range_error("ERROR: Number of symmetry functions set to"
173 "zero, cannot allocate.\n");
174 }
175
176 // Clear all symmetry function related vectors (also for derivatives).
177 G.clear();
178 dEdG.clear();
179 dQdG.clear();
180 dChidG.clear();
181#ifdef N2P2_FULL_SFD_MEMORY
182 dGdxia.clear();
183#endif
184 dGdr.clear();
185
186 // Don't need any allocation of neighbors outside the symmetry function
187 // cutoff.
188 size_t const numSymFuncNeighbors = getStoredMinNumNeighbors(maxCutoffRadius);
189 for (size_t i = 0; i < numSymFuncNeighbors; ++i)
190 {
191 Neighbor& n = neighbors.at(i);
192#ifndef N2P2_NO_SF_CACHE
193 n.cache.clear();
194#endif
195 n.dGdr.clear();
196 }
197
198 // Reset status of symmetry functions and derivatives.
199 hasSymmetryFunctions = false;
201
202 // Resize vectors (derivatives only if requested).
203 G.resize(numSymmetryFunctions, 0.0);
204 if (all)
205 {
206#ifndef N2P2_FULL_SFD_MEMORY
207 if (numSymmetryFunctionDerivatives.size() == 0)
208 {
209 throw range_error("ERROR: Number of symmetry function derivatives"
210 " unset, cannot allocate.\n");
211 }
212#endif
213 if (useChargeNeuron) dEdG.resize(numSymmetryFunctions + 1, 0.0);
214 else dEdG.resize(numSymmetryFunctions, 0.0);
215 dQdG.resize(numSymmetryFunctions, 0.0);
216 dChidG.resize(numSymmetryFunctions,0.0);
217#ifdef N2P2_FULL_SFD_MEMORY
218 dGdxia.resize(numSymmetryFunctions, 0.0);
219#endif
221 }
222 for (size_t i = 0; i < numSymFuncNeighbors; ++i)
223 {
224 Neighbor& n = neighbors.at(i);
225#ifndef N2P2_NO_SF_CACHE
226 n.cache.resize(cacheSizePerElement.at(n.element),
227 -numeric_limits<double>::max());
228#endif
229 if (all)
230 {
231#ifndef N2P2_FULL_SFD_MEMORY
233#else
234 n.dGdr.resize(numSymmetryFunctions);
235#endif
236 }
237 }
238
239 return;
240}
241
242void Atom::free(bool all, double const maxCutoffRadius)
243{
244 size_t const numSymFuncNeighbors = getStoredMinNumNeighbors(maxCutoffRadius);
245 if (all)
246 {
247 G.clear();
248 vector<double>(G).swap(G);
249 hasSymmetryFunctions = false;
250
251#ifndef N2P2_NO_SF_CACHE
252 for (size_t i = 0; i < numSymFuncNeighbors; ++i)
253 {
254 Neighbor& n = neighbors.at(i);
255 n.cache.clear();
256 vector<double>(n.cache).swap(n.cache);
257 }
258#endif
259 }
260
261 dEdG.clear();
262 vector<double>(dEdG).swap(dEdG);
263 dQdG.clear();
264 vector<double>(dQdG).swap(dQdG);
265 dChidG.clear();
266 vector<double>(dChidG).swap(dChidG);
267#ifdef N2P2_FULL_SFD_MEMORY
268 dGdxia.clear();
269 vector<double>(dGdxia).swap(dGdxia);
270#endif
271 dGdr.clear();
272 vector<Vec3D>(dGdr).swap(dGdr);
273
274 for (size_t i = 0; i < numSymFuncNeighbors; ++i)
275 {
276 Neighbor& n = neighbors.at(i);
277 n.dGdr.clear();
278 vector<Vec3D>(n.dGdr).swap(n.dGdr);
279 }
281
282 return;
283}
284
286{
288
289 return;
290}
291
292void Atom::clearNeighborList(size_t const numElements)
293{
294 free(true);
295 numNeighbors = 0;
296 numNeighborsPerElement.resize(numElements, 0);
298 neighborsUnique.clear();
299 vector<size_t>(neighborsUnique).swap(neighborsUnique);
300 neighbors.clear();
301 vector<Atom::Neighbor>(neighbors).swap(neighbors);
302 hasNeighborList = false;
303
304 return;
305}
306
307size_t Atom::calculateNumNeighbors(double const cutoffRadius) const
308{
309 // neighborCutoffs map is only constructed when the neighbor list is sorted.
310 if (unorderedMapContainsKey(neighborCutoffs, cutoffRadius))
311 return neighborCutoffs.at(cutoffRadius);
312 else
313 {
314 size_t numNeighborsLocal = 0;
315
316 for (vector<Neighbor>::const_iterator it = neighbors.begin();
317 it != neighbors.end(); ++it)
318 {
319 if (it->d <= cutoffRadius)
320 {
321 numNeighborsLocal++;
322 }
323 }
324
325 return numNeighborsLocal;
326 }
327}
328
329size_t Atom::getStoredMinNumNeighbors(double const cutoffRadius) const
330{
331 // neighborCutoffs map is only constructed when the neighbor list is sorted.
332 if (unorderedMapContainsKey(neighborCutoffs, cutoffRadius))
333 return neighborCutoffs.at(cutoffRadius);
334 else
335 {
336 return numNeighbors;
337 }
338}
339
340bool Atom::isNeighbor(size_t index) const
341{
342 for (auto const& n : neighbors) if (n.index == index) return true;
343
344 return false;
345}
346
347void Atom::updateError(string const& property,
348 map<string, double>& error,
349 size_t& count) const
350{
351 if (property == "force")
352 {
353 count += 3;
354 error.at("RMSE") += (fRef - f).norm2();
355 error.at("MAE") += (fRef - f).l1norm();
356 }
357 else if (property == "charge")
358 {
359 count++;
360 error.at("RMSE") += (chargeRef - charge) * (chargeRef - charge);
361 error.at("MAE") += fabs(chargeRef - charge);
362 }
363 else
364 {
365 throw runtime_error("ERROR: Unknown property for error update.\n");
366 }
367
368 return;
369}
370
372{
373 Vec3D selfForce{};
374 for (size_t i = 0; i < numSymmetryFunctions; ++i)
375 {
376 selfForce -= dEdG.at(i) * dGdr.at(i);
377 }
378 return selfForce;
379}
380
382 vector<vector<size_t> >
383 const *const tableFull) const
384{
385 Vec3D pairForce{};
386#ifndef N2P2_FULL_SFD_MEMORY
387 if (!tableFull) throw runtime_error(
388 "ERROR: tableFull must not be null pointer");
389 vector<size_t> const& table = tableFull->at(neighbor.element);
390 for (size_t i = 0; i < neighbor.dGdr.size(); ++i)
391 {
392 pairForce -= dEdG.at(table.at(i)) * neighbor.dGdr.at(i);
393 }
394#else
395 for (size_t i = 0; i < numSymmetryFunctions; ++i)
396 {
397 pairForce -= dEdG.at(i) * neighbor.dGdr.at(i);
398 }
399#endif
400 return pairForce;
401}
402
403Vec3D Atom::calculateDChidr(size_t const atomIndexOfR,
404 double const maxCutoffRadius,
405 vector<vector<size_t> > const *const tableFull) const
406{
407#ifndef N2P2_FULL_SFD_MEMORY
408 if (!tableFull) throw runtime_error(
409 "ERROR: tableFull must not be null pointer");
410#endif
411 Vec3D dChidr{};
412 // need to add this case because the loop over the neighbors
413 // does not include the contribution dChi_i/dr_i.
414 if (atomIndexOfR == index)
415 {
416 for (size_t k = 0; k < numSymmetryFunctions; ++k)
417 {
418 dChidr += dChidG.at(k) * dGdr.at(k);
419 }
420 }
421
422 size_t numNeighbors = getStoredMinNumNeighbors(maxCutoffRadius);
423 for (size_t m = 0; m < numNeighbors; ++m)
424 {
425 Atom::Neighbor const& n = neighbors.at(m);
426 // atomIndexOfR must not be larger than ~max(int64_t)/2.
427 if (n.tag == (int64_t)atomIndexOfR)
428 {
429#ifndef N2P2_FULL_SFD_MEMORY
430 vector<size_t> const& table = tableFull->at(n.element);
431 for (size_t k = 0; k < n.dGdr.size(); ++k)
432 {
433 dChidr += dChidG.at(table.at(k)) * n.dGdr.at(k);
434 }
435#else
436 for (size_t k = 0; k < numSymmetryFunctions; ++k)
437 {
438 dChidr += dChidG.at(k) * n.dGdr.at(k);
439 }
440#endif
441 }
442 }
443 return dChidr;
444}
445
446vector<string> Atom::getForcesLines() const
447{
448 vector<string> v;
449 for (size_t i = 0; i < 3; ++i)
450 {
451 v.push_back(strpr("%10zu %10zu %16.8E %16.8E\n",
453 index,
454 fRef[i],
455 f[i]));
456 }
457
458 return v;
459}
460
462{
463 return strpr("%10zu %10zu %16.8E %16.8E\n",
465 index,
466 chargeRef,
467 charge);
468}
469
470vector<string> Atom::info() const
471{
472 vector<string> v;
473
474 v.push_back(strpr("********************************\n"));
475 v.push_back(strpr("ATOM \n"));
476 v.push_back(strpr("********************************\n"));
477 v.push_back(strpr("hasNeighborList : %d\n", hasNeighborList));
478 v.push_back(strpr("hasSymmetryFunctions : %d\n", hasSymmetryFunctions));
479 v.push_back(strpr("hasSymmetryFunctionDerivatives : %d\n", hasSymmetryFunctionDerivatives));
480 v.push_back(strpr("useChargeNeuron : %d\n", useChargeNeuron));
481 v.push_back(strpr("index : %d\n", index));
482 v.push_back(strpr("indexStructure : %d\n", indexStructure));
483 v.push_back(strpr("tag : %" PRId64 "\n", tag));
484 v.push_back(strpr("element : %d\n", element));
485 v.push_back(strpr("numNeighbors : %d\n", numNeighbors));
486 v.push_back(strpr("numNeighborsUnique : %d\n", numNeighborsUnique));
487 v.push_back(strpr("numSymmetryFunctions : %d\n", numSymmetryFunctions));
488 v.push_back(strpr("energy : %16.8E\n", energy));
489 v.push_back(strpr("chi : %16.8E\n", chi));
490 v.push_back(strpr("charge : %16.8E\n", charge));
491 v.push_back(strpr("chargeRef : %16.8E\n", chargeRef));
492 v.push_back(strpr("r : %16.8E %16.8E %16.8E\n", r[0], r[1], r[2]));
493 v.push_back(strpr("f : %16.8E %16.8E %16.8E\n", f[0], f[1], f[2]));
494 v.push_back(strpr("fRef : %16.8E %16.8E %16.8E\n", fRef[0], fRef[1], fRef[2]));
495 v.push_back(strpr("--------------------------------\n"));
496 v.push_back(strpr("neighborsUnique [*] : %d\n", neighborsUnique.size()));
497 v.push_back(strpr("--------------------------------\n"));
498 for (size_t i = 0; i < neighborsUnique.size(); ++i)
499 {
500 v.push_back(strpr("%29d : %d\n", i, neighborsUnique.at(i)));
501 }
502 v.push_back(strpr("--------------------------------\n"));
503 v.push_back(strpr("--------------------------------\n"));
504 v.push_back(strpr("numNeighborsPerElement [*] : %d\n", numNeighborsPerElement.size()));
505 v.push_back(strpr("--------------------------------\n"));
506 for (size_t i = 0; i < numNeighborsPerElement.size(); ++i)
507 {
508 v.push_back(strpr("%29d : %d\n", i, numNeighborsPerElement.at(i)));
509 }
510 v.push_back(strpr("--------------------------------\n"));
511 v.push_back(strpr("--------------------------------\n"));
512 v.push_back(strpr("numSymmetryFunctionDeriv. [*] : %d\n", numSymmetryFunctionDerivatives.size()));
513 v.push_back(strpr("--------------------------------\n"));
514 for (size_t i = 0; i < numSymmetryFunctionDerivatives.size(); ++i)
515 {
516 v.push_back(strpr("%29d : %d\n", i, numSymmetryFunctionDerivatives.at(i)));
517 }
518#ifndef N2P2_NO_SF_CACHE
519 v.push_back(strpr("--------------------------------\n"));
520 v.push_back(strpr("--------------------------------\n"));
521 v.push_back(strpr("cacheSizePerElement [*] : %d\n", cacheSizePerElement.size()));
522 v.push_back(strpr("--------------------------------\n"));
523 for (size_t i = 0; i < cacheSizePerElement.size(); ++i)
524 {
525 v.push_back(strpr("%29d : %d\n", i, cacheSizePerElement.at(i)));
526 }
527#endif
528 v.push_back(strpr("--------------------------------\n"));
529 v.push_back(strpr("--------------------------------\n"));
530 v.push_back(strpr("G [*] : %d\n", G.size()));
531 v.push_back(strpr("--------------------------------\n"));
532 for (size_t i = 0; i < G.size(); ++i)
533 {
534 v.push_back(strpr("%29d : %16.8E\n", i, G.at(i)));
535 }
536 v.push_back(strpr("--------------------------------\n"));
537 v.push_back(strpr("--------------------------------\n"));
538 v.push_back(strpr("dEdG [*] : %d\n", dEdG.size()));
539 v.push_back(strpr("--------------------------------\n"));
540 for (size_t i = 0; i < dEdG.size(); ++i)
541 {
542 v.push_back(strpr("%29d : %16.8E\n", i, dEdG.at(i)));
543 }
544 v.push_back(strpr("--------------------------------\n"));
545 v.push_back(strpr("--------------------------------\n"));
546 v.push_back(strpr("dQdG [*] : %d\n", dQdG.size()));
547 v.push_back(strpr("--------------------------------\n"));
548 for (size_t i = 0; i < dQdG.size(); ++i)
549 {
550 v.push_back(strpr("%29d : %16.8E\n", i, dQdG.at(i)));
551 }
552 v.push_back(strpr("--------------------------------\n"));
553#ifdef N2P2_FULL_SFD_MEMORY
554 v.push_back(strpr("--------------------------------\n"));
555 v.push_back(strpr("dGdxia [*] : %d\n", dGdxia.size()));
556 v.push_back(strpr("--------------------------------\n"));
557 for (size_t i = 0; i < dGdxia.size(); ++i)
558 {
559 v.push_back(strpr("%29d : %16.8E\n", i, dGdxia.at(i)));
560 }
561 v.push_back(strpr("--------------------------------\n"));
562#endif
563 v.push_back(strpr("--------------------------------\n"));
564 v.push_back(strpr("dGdr [*] : %d\n", dGdr.size()));
565 v.push_back(strpr("--------------------------------\n"));
566 for (size_t i = 0; i < dGdr.size(); ++i)
567 {
568 v.push_back(strpr("%29d : %16.8E %16.8E %16.8E\n", i, dGdr.at(i)[0], dGdr.at(i)[1], dGdr.at(i)[2]));
569 }
570 v.push_back(strpr("--------------------------------\n"));
571 v.push_back(strpr("--------------------------------\n"));
572 v.push_back(strpr("neighbors [*] : %d\n", neighbors.size()));
573 v.push_back(strpr("--------------------------------\n"));
574 for (size_t i = 0; i < neighbors.size(); ++i)
575 {
576 v.push_back(strpr("%29d :\n", i));
577 vector<string> vn = neighbors[i].info();
578 v.insert(v.end(), vn.begin(), vn.end());
579 }
580 v.push_back(strpr("--------------------------------\n"));
581 v.push_back(strpr("********************************\n"));
582
583 return v;
584}
585
587 tag (0 ),
588 element (0 ),
589 d (0.0 )
590{
591}
592
594{
595 if (element != rhs.element) return false;
596 if (d != rhs.d ) return false;
597 return true;
598}
599
601{
602 // sorting according to elements deactivated
603 // TODO: resolve this issue for nnp-sfclust
604 // TODO: this breaks pynnp CI tests in test_Neighbor.py
605 //if (element < rhs.element) return true;
606 //else if (element > rhs.element) return false;
607 if (d < rhs.d ) return true;
608 else if (d > rhs.d ) return false;
609 return false;
610}
611
612vector<string> Atom::Neighbor::info() const
613{
614 vector<string> v;
615
616 v.push_back(strpr("********************************\n"));
617 v.push_back(strpr("NEIGHBOR \n"));
618 v.push_back(strpr("********************************\n"));
619 v.push_back(strpr("index : %d\n", index));
620 v.push_back(strpr("tag : %" PRId64 "\n", tag));
621 v.push_back(strpr("element : %d\n", element));
622 v.push_back(strpr("d : %16.8E\n", d));
623 v.push_back(strpr("dr : %16.8E %16.8E %16.8E\n", dr[0], dr[1], dr[2]));
624 v.push_back(strpr("--------------------------------\n"));
625#ifndef N2P2_NO_SF_CACHE
626 v.push_back(strpr("cache [*] : %d\n", cache.size()));
627 v.push_back(strpr("--------------------------------\n"));
628 for (size_t i = 0; i < cache.size(); ++i)
629 {
630 v.push_back(strpr("%29d : %16.8E\n", i, cache.at(i)));
631 }
632 v.push_back(strpr("--------------------------------\n"));
633 v.push_back(strpr("--------------------------------\n"));
634#endif
635 v.push_back(strpr("dGdr [*] : %d\n", dGdr.size()));
636 v.push_back(strpr("--------------------------------\n"));
637 for (size_t i = 0; i < dGdr.size(); ++i)
638 {
639 v.push_back(strpr("%29d : %16.8E %16.8E %16.8E\n", i, dGdr.at(i)[0], dGdr.at(i)[1], dGdr.at(i)[2]));
640 }
641 v.push_back(strpr("--------------------------------\n"));
642 v.push_back(strpr("********************************\n"));
643
644 return v;
645}
Definition Atom.h:29
string strpr(const char *format,...)
String version of printf function.
Definition utility.cpp:90
bool unorderedMapContainsKey(std::unordered_map< T, U > const &stdMap, T key)
Test if unordered map contains specified key.
Definition utility.h:179
Struct to store information on neighbor atoms.
Definition Atom.h:36
Neighbor()
Neighbor constructor, initialize to zero.
Definition Atom.cpp:586
std::size_t index
Index of neighbor atom.
Definition Atom.h:38
std::vector< double > cache
Symmetry function cache (e.g. for cutoffs, compact functions).
Definition Atom.h:49
std::size_t element
Element index of neighbor atom.
Definition Atom.h:42
double d
Distance to neighbor atom.
Definition Atom.h:44
Vec3D dr
Distance vector to neighbor atom.
Definition Atom.h:46
bool operator==(Neighbor const &rhs) const
Overload == operator.
Definition Atom.cpp:593
int64_t tag
Tag of neighbor atom.
Definition Atom.h:40
std::vector< std::string > info() const
Get atom information as a vector of strings.
Definition Atom.cpp:612
std::vector< Vec3D > dGdr
Derivatives of symmetry functions with respect to neighbor coordinates.
Definition Atom.h:60
bool operator<(Neighbor const &rhs) const
Overload < operator.
Definition Atom.cpp:600
std::vector< Neighbor > neighbors
Neighbor array (maximum number defined in macros.h.
Definition Atom.h:170
std::vector< std::string > info() const
Get atom information as a vector of strings.
Definition Atom.cpp:470
std::size_t numSymmetryFunctions
Number of symmetry functions used to describe the atom environment.
Definition Atom.h:113
Vec3D r
Cartesian coordinates.
Definition Atom.h:125
std::vector< double > dEdG
Derivative of atomic energy with respect to symmetry functions.
Definition Atom.h:149
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::vector< double > dQdG
Derivative of atomic charge with respect to symmetry functions.
Definition Atom.h:151
double charge
Atomic charge determined by neural network.
Definition Atom.h:121
bool isNeighbor(std::size_t index) const
Return whether atom is a neighbor.
Definition Atom.cpp:340
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
Vec3D fRef
Reference force vector from data set.
Definition Atom.h:131
bool useChargeNeuron
If an additional charge neuron in the short-range NN is present.
Definition Atom.h:99
std::vector< Vec3D > dGdr
Derivative of symmetry functions with respect to this atom's coordinates.
Definition Atom.h:161
bool NeighborListIsSorted
If the neighbor list is sorted by distance.
Definition Atom.h:93
double chi
Atomic electronegativity determined by neural network.
Definition Atom.h:119
void toPhysicalUnits(double convEnergy, double convLength, double convCharge)
Switch to physical length, energy and charge units.
Definition Atom.cpp:123
void clearNeighborList()
Clear neighbor list.
Definition Atom.cpp:285
void free(bool all, double const maxCutoffRadius=0.0)
Free vectors related to symmetry functions, opposite of allocate().
Definition Atom.cpp:242
std::size_t indexStructure
Index number of structure this atom belongs to.
Definition Atom.h:103
int64_t tag
Tag number of this atom.
Definition Atom.h:105
std::size_t element
Element index of this atom.
Definition Atom.h:107
std::unordered_map< double, size_t > neighborCutoffs
Map stores number of neighbors needed for the corresponding cut-off.
Definition Atom.h:172
void allocate(bool all, double const maxCutoffRadius=0.0)
Allocate vectors related to symmetry functions (G, dEdG).
Definition Atom.cpp:168
void toNormalizedUnits(double convEnergy, double convLength, double convCharge)
Switch to normalized length, energy and charge units.
Definition Atom.cpp:78
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
void updateError(std::string const &property, std::map< std::string, double > &error, std::size_t &count) const
Update property error metrics with data from this atom.
Definition Atom.cpp:347
std::vector< std::size_t > cacheSizePerElement
Cache size for each element.
Definition Atom.h:143
Atom()
Atom constructor, initialize to zero.
Definition Atom.cpp:27
bool hasNeighborList
If the neighbor list has been calculated for this atom.
Definition Atom.h:91
double dEelecdQ
Derivative of electrostatic energy with respect to this atom's charge.
Definition Atom.h:117
double energy
Atomic energy determined by neural network.
Definition Atom.h:115
std::vector< double > G
Symmetry function values.
Definition Atom.h:146
std::vector< std::size_t > neighborsUnique
List of unique neighbor indices (don't count multiple PBC images).
Definition Atom.h:136
double chargeRef
Atomic reference charge.
Definition Atom.h:123
std::vector< std::size_t > numNeighborsPerElement
Number of neighbors per element.
Definition Atom.h:138
std::vector< std::string > getForcesLines() const
Get reference and NN forces for this atoms.
Definition Atom.cpp:446
std::string getChargeLine() const
Get reference and NN charge for this atoms.
Definition Atom.cpp:461
std::size_t numNeighborsUnique
Number of unique neighbor indices (don't count multiple PBC images).
Definition Atom.h:111
std::vector< double > dChidG
Derivative of electronegativity with respect to symmetry functions.
Definition Atom.h:153
std::size_t numNeighbors
Total number of neighbors.
Definition Atom.h:109
Vector in 3 dimensional real space.
Definition Vec3D.h:30