n2p2 - A neural network potential package
nnp::Dataset Class Reference

Collect and process large data sets. More...

#include <Dataset.h>

Inheritance diagram for nnp::Dataset:
Collaboration diagram for nnp::Dataset:

Public Member Functions

 Dataset ()
 Constructor, initialize members. More...
 
 ~Dataset ()
 Destructor. More...
 
void setupMPI ()
 Initialize MPI with MPI_COMM_WORLD. More...
 
void setupMPI (MPI_Comm *communicator)
 Initialize MPI with given communicator. More...
 
void setupRandomNumberGenerator ()
 Initialize random number generator. More...
 
std::size_t getNumStructures (std::ifstream &dataFile)
 Get number of structures in data file. More...
 
int calculateBufferSize (Structure const &structure) const
 Calculate buffer size required to communicate structure via MPI. More...
 
int sendStructure (Structure const &structure, int dest) const
 Send one structure to destination process. More...
 
int recvStructure (Structure *structure, int src)
 Receive one structure from source process. More...
 
int distributeStructures (bool randomize, bool excludeRank0=false, std::string const &fileName="input.data")
 Read data file and distribute structures among processors. More...
 
std::size_t prepareNumericForces (Structure &original, double delta)
 Prepare numeric force check for a single structure. More...
 
void toNormalizedUnits ()
 Switch all structures to normalized units. More...
 
void toPhysicalUnits ()
 Switch all structures to physical units. More...
 
void collectSymmetryFunctionStatistics ()
 Collect symmetry function statistics from all processors. More...
 
void writeSymmetryFunctionScaling (std::string const &fileName="scaling.data")
 Write symmetry function scaling values to file. More...
 
void writeSymmetryFunctionHistograms (std::size_t numBins, std::string fileNameFormat="sf.%03zu.%04zu.histo")
 Calculate and write symmetry function histograms. More...
 
void writeSymmetryFunctionFile (std::string fileName="function.data")
 Write symmetry function legacy file ("function.data"). More...
 
std::size_t writeNeighborHistogram (std::string const &fileNameHisto="neighbors.histo", std::string const &fileNameStructure="neighbors.out")
 Calculate and write neighbor histogram and per-structure statistics. More...
 
void sortNeighborLists ()
 Sort all neighbor lists according to element and distance. More...
 
void writeNeighborLists (std::string const &fileName="neighbor-list.data")
 Write neighbor list file. More...
 
void writeAtomicEnvironmentFile (std::vector< std::vector< std::size_t > > neighCutoff, bool derivatives, std::string const &fileNamePrefix="atomic-env")
 Write atomic environment file. More...
 
void collectError (std::string const &property, std::map< std::string, double > &error, std::size_t &count) const
 Collect error metrics of a property over all MPI procs. More...
 
void combineFiles (std::string filePrefix) const
 Combine individual MPI proc files to one. More...
 
- Public Member Functions inherited from nnp::Mode
 Mode ()
 
void initialize ()
 Write welcome message with version information. More...
 
void loadSettingsFile (std::string const &fileName="input.nn")
 Open settings file and load all keywords into memory. More...
 
void setupGeneric (bool skipNormalize=false)
 Combine multiple setup routines and provide a basic NNP setup. More...
 
void setupNormalization (bool standalone=true)
 Set up normalization. More...
 
virtual void setupElementMap ()
 Set up the element map. More...
 
virtual void setupElements ()
 Set up all Element instances. More...
 
void setupCutoff ()
 Set up cutoff function for all symmetry functions. More...
 
virtual void setupSymmetryFunctions ()
 Set up all symmetry functions. More...
 
void setupSymmetryFunctionScalingNone ()
 Set up "empy" symmetry function scaling. More...
 
virtual void setupSymmetryFunctionScaling (std::string const &fileName="scaling.data")
 Set up symmetry function scaling from file. More...
 
virtual void setupSymmetryFunctionGroups ()
 Set up symmetry function groups. More...
 
virtual void setupSymmetryFunctionCache (bool verbose=false)
 Set up symmetry function cache. More...
 
void setupSymmetryFunctionMemory (bool verbose=false)
 Extract required memory dimensions for symmetry function derivatives. More...
 
void setupSymmetryFunctionStatistics (bool collectStatistics, bool collectExtrapolationWarnings, bool writeExtrapolationWarnings, bool stopOnExtrapolationWarnings)
 Set up symmetry function statistics collection. More...
 
virtual void setupNeuralNetwork ()
 Set up neural networks for all elements. More...
 
virtual void setupNeuralNetworkWeights (std::string const &fileNameFormatShort="weights.%03zu.data", std::string const &fileNameFormatCharge="weightse.%03zu.data")
 Set up neural network weights from files. More...
 
void calculateSymmetryFunctions (Structure &structure, bool const derivatives)
 Calculate all symmetry functions for all atoms in given structure. More...
 
void calculateSymmetryFunctionGroups (Structure &structure, bool const derivatives)
 Calculate all symmetry function groups for all atoms in given structure. More...
 
void calculateAtomicNeuralNetworks (Structure &structure, bool const derivatives)
 Calculate a single atomic neural network for a given atom and nn type. More...
 
void calculateEnergy (Structure &structure) const
 Calculate potential energy for a given structure. More...
 
void calculateCharge (Structure &structure) const
 Calculate total charge for a given structure. More...
 
void calculateForces (Structure &structure) const
 Calculate forces for all atoms in given structure. More...
 
void addEnergyOffset (Structure &structure, bool ref=true)
 Add atomic energy offsets to reference energy. More...
 
void removeEnergyOffset (Structure &structure, bool ref=true)
 Remove atomic energy offsets from reference energy. More...
 
double getEnergyOffset (Structure const &structure) const
 Get atomic energy offset for given structure. More...
 
double getEnergyWithOffset (Structure const &structure, bool ref=true) const
 Add atomic energy offsets and return energy. More...
 
double normalized (std::string const &property, double value) const
 Apply normalization to given property. More...
 
double normalizedEnergy (Structure const &structure, bool ref=true) const
 Apply normalization to given energy of structure. More...
 
double physical (std::string const &property, double value) const
 Undo normalization for a given property. More...
 
double physicalEnergy (Structure const &structure, bool ref=true) const
 Undo normalization for a given energy of structure. More...
 
void convertToNormalizedUnits (Structure &structure) const
 Convert one structure to normalized units. More...
 
void convertToPhysicalUnits (Structure &structure) const
 Convert one structure to physical units. More...
 
std::size_t getNumExtrapolationWarnings () const
 Count total number of extrapolation warnings encountered for all elements and symmetry functions. More...
 
void resetExtrapolationWarnings ()
 Erase all extrapolation warnings and reset counters. More...
 
NNPType getNnpType () const
 Getter for Mode::nnpType. More...
 
double getMeanEnergy () const
 Getter for Mode::meanEnergy. More...
 
double getConvEnergy () const
 Getter for Mode::convEnergy. More...
 
double getConvLength () const
 Getter for Mode::convLength. More...
 
double getMaxCutoffRadius () const
 Getter for Mode::maxCutoffRadius. More...
 
std::size_t getNumElements () const
 Getter for Mode::numElements. More...
 
std::vector< std::size_t > getNumSymmetryFunctions () const
 Get number of symmetry functions per element. More...
 
bool useNormalization () const
 Check if normalization is enabled. More...
 
bool settingsKeywordExists (std::string const &keyword) const
 Check if keyword was found in settings file. More...
 
std::string settingsGetValue (std::string const &keyword) const
 Get value for given keyword in Settings instance. More...
 
std::vector< std::size_t > pruneSymmetryFunctionsRange (double threshold)
 Prune symmetry functions according to their range and write settings file. More...
 
std::vector< std::size_t > pruneSymmetryFunctionsSensitivity (double threshold, std::vector< std::vector< double > > sensitivity)
 Prune symmetry functions with sensitivity analysis data. More...
 
void writePrunedSettingsFile (std::vector< std::size_t > prune, std::string fileName="output.nn") const
 Copy settings file but comment out lines provided. More...
 
void writeSettingsFile (std::ofstream *const &file) const
 Write complete settings file. More...
 

Public Attributes

std::vector< Structurestructures
 All structures in this dataset. More...
 
- Public Attributes inherited from nnp::Mode
ElementMap elementMap
 Global element map, populated by setupElementMap(). More...
 
Log log
 Global log file. More...
 

Protected Attributes

int myRank
 My process ID. More...
 
int numProcs
 Total number of MPI processors. More...
 
std::size_t numStructures
 Total number of structures in dataset. More...
 
std::string myName
 My processor name. More...
 
MPI_Comm comm
 Global MPI communicator. More...
 
gsl_rng * rng
 GSL random number generator (different seed for each MPI process). More...
 
gsl_rng * rngGlobal
 Global GSL random number generator (equal seed for each MPI process). More...
 
- Protected Attributes inherited from nnp::Mode
NNPType nnpType
 
bool normalize
 
bool checkExtrapolationWarnings
 
bool useChargeNN
 
std::size_t numElements
 
std::vector< std::size_t > minNeighbors
 
std::vector< double > minCutoffRadius
 
double maxCutoffRadius
 
double cutoffAlpha
 
double meanEnergy
 
double convEnergy
 
double convLength
 
Settings settings
 
SymFnc::ScalingType scalingType
 
CutoffFunction::CutoffType cutoffType
 
std::vector< Elementelements
 

Additional Inherited Members

- Public Types inherited from nnp::Mode
enum class  NNPType { SHORT_ONLY , SHORT_CHARGE_NN }
 
- Protected Member Functions inherited from nnp::Mode
void readNeuralNetworkWeights (std::string const &type, std::string const &fileName)
 Read in weights for a specific type of neural network. More...
 

Detailed Description

Collect and process large data sets.

Definition at line 34 of file Dataset.h.

Constructor & Destructor Documentation

◆ Dataset()

Dataset::Dataset ( )

Constructor, initialize members.

Definition at line 36 of file Dataset.cpp.

36 : Mode(),
37 myRank (0 ),
38 numProcs (0 ),
39 numStructures(0 ),
40 myName ("" ),
41 rng (NULL),
42 rngGlobal (NULL)
43{
44}
std::size_t numStructures
Total number of structures in dataset.
Definition: Dataset.h:203
gsl_rng * rng
GSL random number generator (different seed for each MPI process).
Definition: Dataset.h:209
int numProcs
Total number of MPI processors.
Definition: Dataset.h:201
std::string myName
My processor name.
Definition: Dataset.h:205
int myRank
My process ID.
Definition: Dataset.h:199
gsl_rng * rngGlobal
Global GSL random number generator (equal seed for each MPI process).
Definition: Dataset.h:211
Mode()
Definition: Mode.cpp:37

◆ ~Dataset()

Dataset::~Dataset ( )

Destructor.

Definition at line 46 of file Dataset.cpp.

47{
48 if (rng != NULL) gsl_rng_free(rng);
49 if (rngGlobal != NULL) gsl_rng_free(rngGlobal);
50}

References rng, and rngGlobal.

Member Function Documentation

◆ setupMPI() [1/2]

void Dataset::setupMPI ( )

Initialize MPI with MPI_COMM_WORLD.

Definition at line 52 of file Dataset.cpp.

53{
54 MPI_Comm tmpComm;
55 MPI_Comm_dup(MPI_COMM_WORLD, &tmpComm);
56 setupMPI(&tmpComm);
57 MPI_Comm_free(&tmpComm);
58
59 return;
60}
void setupMPI()
Initialize MPI with MPI_COMM_WORLD.
Definition: Dataset.cpp:52

References setupMPI().

Referenced by main(), and setupMPI().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setupMPI() [2/2]

void Dataset::setupMPI ( MPI_Comm *  communicator)

Initialize MPI with given communicator.

Parameters
[in]communicatorProvided communicator which should be used.

Definition at line 62 of file Dataset.cpp.

63{
64 log << "\n";
65 log << "*** SETUP: MPI **************************"
66 "**************************************\n";
67 log << "\n";
68
69 int bufferSize = 0;
70 char line[MPI_MAX_PROCESSOR_NAME];
71 MPI_Status ms;
72
73 MPI_Comm_dup(*communicator, &comm);
74 MPI_Comm_rank(comm, &myRank);
75 MPI_Comm_size(comm, &numProcs);
76 MPI_Get_processor_name(line, &bufferSize);
77 myName = line;
78
79 log << strpr("Number of processors: %d\n", numProcs);
80 log << strpr("Process %d of %d (rank %d): %s\n",
81 myRank + 1,
83 myRank,
84 myName.c_str());
85 for (int i = 1; i < numProcs; ++i)
86 {
87 if (myRank == 0)
88 {
89 MPI_Recv(&bufferSize, 1, MPI_INT, i, 0, comm, &ms);
90 MPI_Recv(line, bufferSize, MPI_CHAR, i, 0, comm, &ms);
91 log << strpr("Process %d of %d (rank %d): %s\n",
92 i + 1,
94 i,
95 line);
96 }
97 else if (myRank == i)
98 {
99 MPI_Send(&bufferSize, 1, MPI_INT, 0, 0, comm);
100 MPI_Send(line, bufferSize, MPI_CHAR, 0, 0, comm);
101 }
102 }
103
104 log << "*****************************************"
105 "**************************************\n";
106
107 return;
108}
MPI_Comm comm
Global MPI communicator.
Definition: Dataset.h:207
Log log
Global log file.
Definition: Mode.h:510
string strpr(const char *format,...)
String version of printf function.
Definition: utility.cpp:90

References comm, nnp::Mode::log, myName, myRank, numProcs, and nnp::strpr().

Here is the call graph for this function:

◆ setupRandomNumberGenerator()

void Dataset::setupRandomNumberGenerator ( )

Initialize random number generator.

CAUTION: MPI communicator required, call setupMPI() before!

Definition at line 110 of file Dataset.cpp.

111{
112 log << "\n";
113 log << "*** SETUP: RANDOM NUMBER GENERATOR ******"
114 "**************************************\n";
115 log << "\n";
116
117 // Get random seed from settings file.
118 unsigned long seed = atoi(settings["random_seed"].c_str());
119 unsigned long seedGlobal = 0;
120
121 if (myRank == 0)
122 {
123 log << strpr("Random number generator seed: %d\n", seed);
124 if (seed == 0)
125 {
126 log << "WARNING: Seed set to 0. This is a special value for the "
127 "gsl_rng_set() routine (see GSL docs).\n";
128 }
129 // Initialize personal RNG of process 0 with the given seed.
130 rng = gsl_rng_alloc(gsl_rng_mt19937);
131 gsl_rng_set(rng, seed);
132 log << strpr("Seed for rank %d: %lu\n", 0, seed);
133 for (int i = 1; i < numProcs; ++i)
134 {
135 // Get new seeds for all remaining processes with the RNG.
136 seed = gsl_rng_get(rng);
137 log << strpr("Seed for rank %d: %lu\n", i, seed);
138 MPI_Send(&seed, 1, MPI_UNSIGNED_LONG, i, 0, comm);
139 }
140 // Set seed for global RNG.
141 seedGlobal = gsl_rng_get(rng);
142 }
143 else
144 {
145 // Receive seed for personal RNG.
146 MPI_Status ms;
147 MPI_Recv(&seed, 1, MPI_UNSIGNED_LONG, 0, 0, comm, &ms);
148 log << strpr("Seed for rank %d: %lu\n", myRank, seed);
149 rng = gsl_rng_alloc(gsl_rng_taus);
150 gsl_rng_set(rng, seed);
151 }
152 // Rank 0 broadcasts global seed.
153 MPI_Bcast(&seedGlobal, 1, MPI_UNSIGNED_LONG, 0, comm);
154 log << strpr("Seed for global RNG: %lu\n", seedGlobal);
155 // All processes initialize global RNG.
156 rngGlobal = gsl_rng_alloc(gsl_rng_taus);
157 gsl_rng_set(rngGlobal, seedGlobal);
158
159 log << "*****************************************"
160 "**************************************\n";
161
162 return;
163}
Settings settings
Definition: Mode.h:525

References comm, nnp::Mode::log, myRank, numProcs, rng, rngGlobal, nnp::Mode::settings, and nnp::strpr().

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getNumStructures()

size_t Dataset::getNumStructures ( std::ifstream &  dataFile)

Get number of structures in data file.

Parameters
[in,out]dataFileData file name.
Returns
Number of structures.

Definition at line 707 of file Dataset.cpp.

708{
709 size_t count = 0;
710 string line;
711 vector<string> splitLine;
712
713 while (getline(dataFile, line))
714 {
715 splitLine = split(reduce(line));
716 if (splitLine.at(0) == "begin") count++;
717 }
718
719 return count;
720}
vector< string > split(string const &input, char delimiter)
Split string at each delimiter.
Definition: utility.cpp:33
string reduce(string const &line, string const &whitespace, string const &fill)
Replace multiple whitespaces with fill.
Definition: utility.cpp:60

References nnp::reduce(), and nnp::split().

Referenced by distributeStructures(), and main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculateBufferSize()

int Dataset::calculateBufferSize ( Structure const &  structure) const

Calculate buffer size required to communicate structure via MPI.

Parameters
[in]structureInput structure.
Returns
Buffer size in bytes.

Definition at line 165 of file Dataset.cpp.

166{
167 int bs = 0; // Send buffer size.
168 int is = 0; // int size.
169 int i64s = 0; // int64_t size.
170 int ss = 0; // size_t size.
171 int ds = 0; // double size.
172 int cs = 0; // char size.
173 Structure const& s = structure; // Shortcut for structure.
174
175 MPI_Pack_size(1, MPI_INT , comm, &is );
176 MPI_Pack_size(1, MPI_INT64_T, comm, &i64s);
177 MPI_Pack_size(1, MPI_SIZE_T , comm, &ss );
178 MPI_Pack_size(1, MPI_DOUBLE , comm, &ds );
179 MPI_Pack_size(1, MPI_CHAR , comm, &cs );
180
181 // Structure
182 bs += 5 * cs + 4 * ss + 4 * is + 5 * ds;
183 // Structure.comment
184 bs += ss;
185 bs += (s.comment.length() + 1) * cs;
186 // Structure.box
187 bs += 9 * ds;
188 // Structure.invbox
189 bs += 9 * ds;
190 // Structure.numAtomsPerElement
191 bs += ss;
192 bs += s.numAtomsPerElement.size() * ss;
193 // Structure.atoms
194 bs += ss;
195 bs += s.atoms.size() * (4 * cs + 6 * ss + i64s + 3 * ds + 3 * 3 * ds);
196 for (vector<Atom>::const_iterator it = s.atoms.begin();
197 it != s.atoms.end(); ++it)
198 {
199 // Atom.neighborsUnique
200 bs += ss;
201 bs += it->neighborsUnique.size() * ss;
202 // Atom.numNeighborsPerElement
203 bs += ss;
204 bs += it->numNeighborsPerElement.size() * ss;
205 // Atom.numSymmetryFunctionDerivatives
206 bs += ss;
207 bs += it->numSymmetryFunctionDerivatives.size() * ss;
208#ifndef N2P2_NO_SF_CACHE
209 // Atom.cacheSizePerElement
210 bs += ss;
211 bs += it->cacheSizePerElement.size() * ss;
212#endif
213 // Atom.G
214 bs += ss;
215 bs += it->G.size() * ds;
216 // Atom.dEdG
217 bs += ss;
218 bs += it->dEdG.size() * ds;
219 // Atom.dQdG
220 bs += ss;
221 bs += it->dQdG.size() * ds;
222#ifdef N2P2_FULL_SFD_MEMORY
223 // Atom.dGdxia
224 bs += ss;
225 bs += it->dGdxia.size() * ds;
226#endif
227 // Atom.dGdr
228 bs += ss;
229 bs += it->dGdr.size() * 3 * ds;
230 // Atom.neighbors
231 bs += ss;
232 for (vector<Atom::Neighbor>::const_iterator it2 =
233 it->neighbors.begin(); it2 != it->neighbors.end(); ++it2)
234 {
235 // Neighbor
236 bs += 2 * ss + i64s + ds + 3 * ds;
237#ifndef N2P2_NO_SF_CACHE
238 // Neighbor.cache
239 bs += ss;
240 bs += it2->cache.size() * ds;
241#endif
242 // Neighbor.dGdr
243 bs += ss;
244 bs += it2->dGdr.size() * 3 * ds;
245 }
246 }
247
248 return bs;
249}
#define MPI_SIZE_T
Definition: mpi-extra.h:22
Storage for one atomic configuration.
Definition: Structure.h:34
std::vector< std::size_t > numAtomsPerElement
Number of atoms of each element in this structure.
Definition: Structure.h:94
std::string comment
Structure comment.
Definition: Structure.h:88
std::vector< Atom > atoms
Vector of all atoms in this structure.
Definition: Structure.h:96

References nnp::Structure::atoms, comm, nnp::Structure::comment, MPI_SIZE_T, and nnp::Structure::numAtomsPerElement.

Referenced by main(), and sendStructure().

Here is the caller graph for this function:

◆ sendStructure()

int Dataset::sendStructure ( Structure const &  structure,
int  dest 
) const

Send one structure to destination process.

Parameters
[in]structureSource structure.
[in]destMPI rank of destination process.
Returns
Number of bytes sent.

Definition at line 251 of file Dataset.cpp.

252{
253 unsigned char* buf = 0; // Send buffer.
254 int bs = 0; // Send buffer size.
255 int p = 0; // Send buffer position.
256 int ts = 0; // Size for temporary stuff.
257 Structure const& s = structure; // Shortcut for structure.
258
259 bs = calculateBufferSize(s);
260 buf = new unsigned char[bs];
261
262 // Structure
263 MPI_Pack(&(s.isPeriodic ), 1, MPI_CHAR , buf, bs, &p, comm);
264 MPI_Pack(&(s.isTriclinic ), 1, MPI_CHAR , buf, bs, &p, comm);
265 MPI_Pack(&(s.hasNeighborList ), 1, MPI_CHAR , buf, bs, &p, comm);
266 MPI_Pack(&(s.hasSymmetryFunctions ), 1, MPI_CHAR , buf, bs, &p, comm);
267 MPI_Pack(&(s.hasSymmetryFunctionDerivatives), 1, MPI_CHAR , buf, bs, &p, comm);
268 MPI_Pack(&(s.index ), 1, MPI_SIZE_T, buf, bs, &p, comm);
269 MPI_Pack(&(s.numAtoms ), 1, MPI_SIZE_T, buf, bs, &p, comm);
270 MPI_Pack(&(s.numElements ), 1, MPI_SIZE_T, buf, bs, &p, comm);
271 MPI_Pack(&(s.numElementsPresent ), 1, MPI_SIZE_T, buf, bs, &p, comm);
272 MPI_Pack(&(s.pbc ), 3, MPI_INT , buf, bs, &p, comm);
273 MPI_Pack(&(s.energy ), 1, MPI_DOUBLE, buf, bs, &p, comm);
274 MPI_Pack(&(s.energyRef ), 1, MPI_DOUBLE, buf, bs, &p, comm);
275 MPI_Pack(&(s.charge ), 1, MPI_DOUBLE, buf, bs, &p, comm);
276 MPI_Pack(&(s.chargeRef ), 1, MPI_DOUBLE, buf, bs, &p, comm);
277 MPI_Pack(&(s.volume ), 1, MPI_DOUBLE, buf, bs, &p, comm);
278 MPI_Pack(&(s.sampleType ), 1, MPI_INT , buf, bs, &p, comm);
279
280 // Strucuture.comment
281 ts = s.comment.length() + 1;
282 MPI_Pack(&ts, 1, MPI_SIZE_T, buf, bs, &p, comm);
283 MPI_Pack(s.comment.c_str(), ts, MPI_CHAR, buf, bs, &p, comm);
284
285 // Structure.box
286 for (size_t i = 0; i < 3; ++i)
287 {
288 MPI_Pack(s.box[i].r, 3, MPI_DOUBLE, buf, bs, &p, comm);
289 }
290
291 // Structure.invbox
292 for (size_t i = 0; i < 3; ++i)
293 {
294 MPI_Pack(s.invbox[i].r, 3, MPI_DOUBLE, buf, bs, &p, comm);
295 }
296
297 // Structure.numAtomsPerElement
298 ts = s.numAtomsPerElement.size();
299 MPI_Pack(&ts, 1, MPI_SIZE_T, buf, bs, &p, comm);
300 if (ts > 0)
301 {
302 MPI_Pack(&(s.numAtomsPerElement.front()), ts, MPI_SIZE_T, buf, bs, &p, comm);
303 }
304
305 // Structure.atoms
306 ts = s.atoms.size();
307 MPI_Pack(&ts, 1, MPI_SIZE_T, buf, bs, &p, comm);
308 if (ts > 0)
309 {
310 for (vector<Atom>::const_iterator it = s.atoms.begin();
311 it != s.atoms.end(); ++it)
312 {
313 // Atom
314 MPI_Pack(&(it->hasNeighborList ), 1, MPI_CHAR , buf, bs, &p, comm);
315 MPI_Pack(&(it->hasSymmetryFunctions ), 1, MPI_CHAR , buf, bs, &p, comm);
316 MPI_Pack(&(it->hasSymmetryFunctionDerivatives), 1, MPI_CHAR , buf, bs, &p, comm);
317 MPI_Pack(&(it->useChargeNeuron ), 1, MPI_CHAR , buf, bs, &p, comm);
318 MPI_Pack(&(it->index ), 1, MPI_SIZE_T , buf, bs, &p, comm);
319 MPI_Pack(&(it->indexStructure ), 1, MPI_SIZE_T , buf, bs, &p, comm);
320 MPI_Pack(&(it->tag ), 1, MPI_INT64_T, buf, bs, &p, comm);
321 MPI_Pack(&(it->element ), 1, MPI_SIZE_T , buf, bs, &p, comm);
322 MPI_Pack(&(it->numNeighbors ), 1, MPI_SIZE_T , buf, bs, &p, comm);
323 MPI_Pack(&(it->numNeighborsUnique ), 1, MPI_SIZE_T , buf, bs, &p, comm);
324 MPI_Pack(&(it->numSymmetryFunctions ), 1, MPI_SIZE_T , buf, bs, &p, comm);
325 MPI_Pack(&(it->energy ), 1, MPI_DOUBLE , buf, bs, &p, comm);
326 MPI_Pack(&(it->charge ), 1, MPI_DOUBLE , buf, bs, &p, comm);
327 MPI_Pack(&(it->chargeRef ), 1, MPI_DOUBLE , buf, bs, &p, comm);
328 MPI_Pack(&(it->r.r ), 3, MPI_DOUBLE , buf, bs, &p, comm);
329 MPI_Pack(&(it->f.r ), 3, MPI_DOUBLE , buf, bs, &p, comm);
330 MPI_Pack(&(it->fRef.r ), 3, MPI_DOUBLE , buf, bs, &p, comm);
331
332 // Atom.neighborsUnique
333 size_t ts2 = it->neighborsUnique.size();
334 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
335 if (ts2 > 0)
336 {
337 MPI_Pack(&(it->neighborsUnique.front()), ts2, MPI_SIZE_T, buf, bs, &p, comm);
338 }
339
340 // Atom.numNeighborsPerElement
341 ts2 = it->numNeighborsPerElement.size();
342 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
343 if (ts2 > 0)
344 {
345 MPI_Pack(&(it->numNeighborsPerElement.front()), ts2, MPI_SIZE_T, buf, bs, &p, comm);
346 }
347
348 // Atom.numSymmetryFunctionDerivatives
349 ts2 = it->numSymmetryFunctionDerivatives.size();
350 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
351 if (ts2 > 0)
352 {
353 MPI_Pack(&(it->numSymmetryFunctionDerivatives.front()), ts2, MPI_SIZE_T, buf, bs, &p, comm);
354 }
355
356#ifndef N2P2_NO_SF_CACHE
357 // Atom.cacheSizePerElement
358 ts2 = it->cacheSizePerElement.size();
359 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
360 if (ts2 > 0)
361 {
362 MPI_Pack(&(it->cacheSizePerElement.front()), ts2, MPI_SIZE_T, buf, bs, &p, comm);
363 }
364#endif
365
366 // Atom.G
367 ts2 = it->G.size();
368 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
369 if (ts2 > 0)
370 {
371 MPI_Pack(&(it->G.front()), ts2, MPI_DOUBLE, buf, bs, &p, comm);
372 }
373
374 // Atom.dEdG
375 ts2 = it->dEdG.size();
376 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
377 if (ts2 > 0)
378 {
379 MPI_Pack(&(it->dEdG.front()), ts2, MPI_DOUBLE, buf, bs, &p, comm);
380 }
381
382 // Atom.dQdG
383 ts2 = it->dQdG.size();
384 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
385 if (ts2 > 0)
386 {
387 MPI_Pack(&(it->dQdG.front()), ts2, MPI_DOUBLE, buf, bs, &p, comm);
388 }
389
390#ifdef N2P2_FULL_SFD_MEMORY
391 // Atom.dGdxia
392 ts2 = it->dGdxia.size();
393 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
394 if (ts2 > 0)
395 {
396 MPI_Pack(&(it->dGdxia.front()), ts2, MPI_DOUBLE, buf, bs, &p, comm);
397 }
398#endif
399
400 // Atom.dGdr
401 ts2 = it->dGdr.size();
402 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
403 if (ts2 > 0)
404 {
405 for (vector<Vec3D>::const_iterator it2 = it->dGdr.begin();
406 it2 != it->dGdr.end(); ++it2)
407 {
408 MPI_Pack(it2->r, 3, MPI_DOUBLE, buf, bs, &p, comm);
409 }
410 }
411
412 // Atom.neighbors
413 ts2 = it->neighbors.size();
414 MPI_Pack(&ts2, 1, MPI_SIZE_T, buf, bs, &p, comm);
415 if (ts2 > 0)
416 {
417 for (vector<Atom::Neighbor>::const_iterator it2 =
418 it->neighbors.begin(); it2 != it->neighbors.end(); ++it2)
419 {
420 // Neighbor
421 MPI_Pack(&(it2->index ), 1, MPI_SIZE_T , buf, bs, &p, comm);
422 MPI_Pack(&(it2->tag ), 1, MPI_INT64_T, buf, bs, &p, comm);
423 MPI_Pack(&(it2->element ), 1, MPI_SIZE_T , buf, bs, &p, comm);
424 MPI_Pack(&(it2->d ), 1, MPI_DOUBLE , buf, bs, &p, comm);
425 MPI_Pack( it2->dr.r , 3, MPI_DOUBLE , buf, bs, &p, comm);
426
427 size_t ts3 = 0;
428#ifndef N2P2_NO_SF_CACHE
429 // Neighbor.cache
430 ts3 = it2->cache.size();
431 MPI_Pack(&ts3, 1, MPI_SIZE_T, buf, bs, &p, comm);
432 if (ts3 > 0)
433 {
434 MPI_Pack(&(it2->cache.front()), ts3, MPI_DOUBLE, buf, bs, &p, comm);
435 }
436#endif
437
438 // Neighbor.dGdr
439 ts3 = it2->dGdr.size();
440 MPI_Pack(&ts3, 1, MPI_SIZE_T, buf, bs, &p, comm);
441 if (ts3 > 0)
442 {
443 for (vector<Vec3D>::const_iterator it3 =
444 it2->dGdr.begin(); it3 != it2->dGdr.end(); ++it3)
445 {
446 MPI_Pack(it3->r, 3, MPI_DOUBLE, buf, bs, &p, comm);
447 }
448 }
449 }
450 }
451 }
452 }
453
454 MPI_Send(&bs, 1, MPI_INT, dest, 0, comm);
455 MPI_Send(buf, bs, MPI_PACKED, dest, 0, comm);
456
457 delete[] buf;
458
459 return bs;
460}
int calculateBufferSize(Structure const &structure) const
Calculate buffer size required to communicate structure via MPI.
Definition: Dataset.cpp:165
size_t p
Definition: nnp-cutoff.cpp:33
Vec3D invbox[3]
Inverse simulation box vectors.
Definition: Structure.h:92
Vec3D box[3]
Simulation box vectors.
Definition: Structure.h:90
bool isTriclinic
If the simulation box is triclinic.
Definition: Structure.h:58
bool isPeriodic
If periodic boundary conditions apply.
Definition: Structure.h:56
double charge
Charge determined by neural network potential.
Definition: Structure.h:80
std::size_t index
Index number of this structure.
Definition: Structure.h:66
double chargeRef
Reference charge.
Definition: Structure.h:82
SampleType sampleType
Sample type (training or test set).
Definition: Structure.h:86
bool hasSymmetryFunctionDerivatives
If symmetry function derivatives are saved for each atom.
Definition: Structure.h:64
double energy
Potential energy determined by neural network.
Definition: Structure.h:76
double energyRef
Reference potential energy.
Definition: Structure.h:78
int pbc[3]
Number of PBC images necessary in each direction.
Definition: Structure.h:74
double volume
Simulation box volume.
Definition: Structure.h:84
std::size_t numAtoms
Total number of atoms present in this structure.
Definition: Structure.h:68
std::size_t numElements
Global number of elements (all structures).
Definition: Structure.h:70
std::size_t numElementsPresent
Number of elements present in this structure.
Definition: Structure.h:72
bool hasNeighborList
If the neighbor list has been calculated.
Definition: Structure.h:60
bool hasSymmetryFunctions
If symmetry function values are saved for each atom.
Definition: Structure.h:62
double r[3]
cartesian coordinates.
Definition: Vec3D.h:31

References nnp::Structure::atoms, nnp::Structure::box, calculateBufferSize(), nnp::Structure::charge, nnp::Structure::chargeRef, comm, nnp::Structure::comment, nnp::Structure::energy, nnp::Structure::energyRef, nnp::Structure::hasNeighborList, nnp::Structure::hasSymmetryFunctionDerivatives, nnp::Structure::hasSymmetryFunctions, nnp::Structure::index, nnp::Structure::invbox, nnp::Structure::isPeriodic, nnp::Structure::isTriclinic, MPI_SIZE_T, nnp::Structure::numAtoms, nnp::Structure::numAtomsPerElement, nnp::Structure::numElements, nnp::Structure::numElementsPresent, p, nnp::Structure::pbc, nnp::Vec3D::r, nnp::Structure::sampleType, and nnp::Structure::volume.

Referenced by distributeStructures(), and prepareNumericForces().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ recvStructure()

int Dataset::recvStructure ( Structure structure,
int  src 
)

Receive one structure from source process.

Parameters
[in]structureSource structure.
[in]srcMPI rank of source process.
Returns
Number of bytes received.

Definition at line 462 of file Dataset.cpp.

463{
464 unsigned char* buf = 0; // Receive buffer.
465 int bs = 0; // Receive buffer size.
466 int p = 0; // Receive buffer position.
467 int ts = 0; // Size for temporary stuff.
468 Structure* const s = structure; // Shortcut for structure.
469 MPI_Status ms;
470
471 // Receive buffer size and extract source.
472 MPI_Recv(&bs, 1, MPI_INT, src, 0, comm, &ms);
473 src = ms.MPI_SOURCE;
474
475 buf = new unsigned char[bs];
476
477 MPI_Recv(buf, bs, MPI_PACKED, src, 0, comm, &ms);
478
479 // Structure
480 MPI_Unpack(buf, bs, &p, &(s->isPeriodic ), 1, MPI_CHAR , comm);
481 MPI_Unpack(buf, bs, &p, &(s->isTriclinic ), 1, MPI_CHAR , comm);
482 MPI_Unpack(buf, bs, &p, &(s->hasNeighborList ), 1, MPI_CHAR , comm);
483 MPI_Unpack(buf, bs, &p, &(s->hasSymmetryFunctions ), 1, MPI_CHAR , comm);
484 MPI_Unpack(buf, bs, &p, &(s->hasSymmetryFunctionDerivatives), 1, MPI_CHAR , comm);
485 MPI_Unpack(buf, bs, &p, &(s->index ), 1, MPI_SIZE_T, comm);
486 MPI_Unpack(buf, bs, &p, &(s->numAtoms ), 1, MPI_SIZE_T, comm);
487 MPI_Unpack(buf, bs, &p, &(s->numElements ), 1, MPI_SIZE_T, comm);
488 MPI_Unpack(buf, bs, &p, &(s->numElementsPresent ), 1, MPI_SIZE_T, comm);
489 MPI_Unpack(buf, bs, &p, &(s->pbc ), 3, MPI_INT , comm);
490 MPI_Unpack(buf, bs, &p, &(s->energy ), 1, MPI_DOUBLE, comm);
491 MPI_Unpack(buf, bs, &p, &(s->energyRef ), 1, MPI_DOUBLE, comm);
492 MPI_Unpack(buf, bs, &p, &(s->charge ), 1, MPI_DOUBLE, comm);
493 MPI_Unpack(buf, bs, &p, &(s->chargeRef ), 1, MPI_DOUBLE, comm);
494 MPI_Unpack(buf, bs, &p, &(s->volume ), 1, MPI_DOUBLE, comm);
495 MPI_Unpack(buf, bs, &p, &(s->sampleType ), 1, MPI_INT , comm);
496
497 // Strucuture.comment
498 ts = 0;
499 MPI_Unpack(buf, bs, &p, &ts, 1, MPI_SIZE_T, comm);
500 char* comment = new char[ts];
501 MPI_Unpack(buf, bs, &p, comment, ts, MPI_CHAR, comm);
502 s->comment = comment;
503 delete[] comment;
504
505 // Structure.box
506 for (size_t i = 0; i < 3; ++i)
507 {
508 MPI_Unpack(buf, bs, &p, s->box[i].r, 3, MPI_DOUBLE, comm);
509 }
510
511 // Structure.invbox
512 for (size_t i = 0; i < 3; ++i)
513 {
514 MPI_Unpack(buf, bs, &p, s->invbox[i].r, 3, MPI_DOUBLE, comm);
515 }
516
517 // Structure.numAtomsPerElement
518 ts = 0;
519 MPI_Unpack(buf, bs, &p, &ts, 1, MPI_SIZE_T, comm);
520 if (ts > 0)
521 {
522 s->numAtomsPerElement.clear();
523 s->numAtomsPerElement.resize(ts, 0);
524 MPI_Unpack(buf, bs, &p, &(s->numAtomsPerElement.front()), ts, MPI_SIZE_T, comm);
525 }
526
527 // Structure.atoms
528 ts = 0;
529 MPI_Unpack(buf, bs, &p, &ts, 1, MPI_SIZE_T, comm);
530 if (ts > 0)
531 {
532 s->atoms.clear();
533 s->atoms.resize(ts);
534 for (vector<Atom>::iterator it = s->atoms.begin();
535 it != s->atoms.end(); ++it)
536 {
537 // Atom
538 MPI_Unpack(buf, bs, &p, &(it->hasNeighborList ), 1, MPI_CHAR , comm);
539 MPI_Unpack(buf, bs, &p, &(it->hasSymmetryFunctions ), 1, MPI_CHAR , comm);
540 MPI_Unpack(buf, bs, &p, &(it->hasSymmetryFunctionDerivatives), 1, MPI_CHAR , comm);
541 MPI_Unpack(buf, bs, &p, &(it->useChargeNeuron ), 1, MPI_CHAR , comm);
542 MPI_Unpack(buf, bs, &p, &(it->index ), 1, MPI_SIZE_T , comm);
543 MPI_Unpack(buf, bs, &p, &(it->indexStructure ), 1, MPI_SIZE_T , comm);
544 MPI_Unpack(buf, bs, &p, &(it->tag ), 1, MPI_INT64_T, comm);
545 MPI_Unpack(buf, bs, &p, &(it->element ), 1, MPI_SIZE_T , comm);
546 MPI_Unpack(buf, bs, &p, &(it->numNeighbors ), 1, MPI_SIZE_T , comm);
547 MPI_Unpack(buf, bs, &p, &(it->numNeighborsUnique ), 1, MPI_SIZE_T , comm);
548 MPI_Unpack(buf, bs, &p, &(it->numSymmetryFunctions ), 1, MPI_SIZE_T , comm);
549 MPI_Unpack(buf, bs, &p, &(it->energy ), 1, MPI_DOUBLE , comm);
550 MPI_Unpack(buf, bs, &p, &(it->charge ), 1, MPI_DOUBLE , comm);
551 MPI_Unpack(buf, bs, &p, &(it->chargeRef ), 1, MPI_DOUBLE , comm);
552 MPI_Unpack(buf, bs, &p, &(it->r.r ), 3, MPI_DOUBLE , comm);
553 MPI_Unpack(buf, bs, &p, &(it->f.r ), 3, MPI_DOUBLE , comm);
554 MPI_Unpack(buf, bs, &p, &(it->fRef.r ), 3, MPI_DOUBLE , comm);
555
556 // Atom.neighborsUnique
557 size_t ts2 = 0;
558 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
559 if (ts2 > 0)
560 {
561 it->neighborsUnique.clear();
562 it->neighborsUnique.resize(ts2, 0);
563 MPI_Unpack(buf, bs, &p, &(it->neighborsUnique.front()), ts2, MPI_SIZE_T, comm);
564 }
565
566 // Atom.numNeighborsPerElement
567 ts2 = 0;
568 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
569 if (ts2 > 0)
570 {
571 it->numNeighborsPerElement.clear();
572 it->numNeighborsPerElement.resize(ts2, 0);
573 MPI_Unpack(buf, bs, &p, &(it->numNeighborsPerElement.front()), ts2, MPI_SIZE_T, comm);
574 }
575
576 // Atom.numSymmetryFunctionDerivatives
577 ts2 = 0;
578 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
579 if (ts2 > 0)
580 {
581 it->numSymmetryFunctionDerivatives.clear();
582 it->numSymmetryFunctionDerivatives.resize(ts2, 0);
583 MPI_Unpack(buf, bs, &p, &(it->numSymmetryFunctionDerivatives.front()), ts2, MPI_SIZE_T, comm);
584 }
585
586#ifndef N2P2_NO_SF_CACHE
587 // Atom.cacheSizePerElement
588 ts2 = 0;
589 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
590 if (ts2 > 0)
591 {
592 it->cacheSizePerElement.clear();
593 it->cacheSizePerElement.resize(ts2, 0);
594 MPI_Unpack(buf, bs, &p, &(it->cacheSizePerElement.front()), ts2, MPI_SIZE_T, comm);
595 }
596#endif
597
598 // Atom.G
599 ts2 = 0;
600 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
601 if (ts2 > 0)
602 {
603 it->G.clear();
604 it->G.resize(ts2, 0.0);
605 MPI_Unpack(buf, bs, &p, &(it->G.front()), ts2, MPI_DOUBLE, comm);
606 }
607
608 // Atom.dEdG
609 ts2 = 0;
610 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
611 if (ts2 > 0)
612 {
613 it->dEdG.clear();
614 it->dEdG.resize(ts2, 0.0);
615 MPI_Unpack(buf, bs, &p, &(it->dEdG.front()), ts2, MPI_DOUBLE, comm);
616 }
617
618 // Atom.dQdG
619 ts2 = 0;
620 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
621 if (ts2 > 0)
622 {
623 it->dQdG.clear();
624 it->dQdG.resize(ts2, 0.0);
625 MPI_Unpack(buf, bs, &p, &(it->dQdG.front()), ts2, MPI_DOUBLE, comm);
626 }
627
628#ifdef N2P2_FULL_SFD_MEMORY
629 // Atom.dGdxia
630 ts2 = 0;
631 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
632 if (ts2 > 0)
633 {
634 it->dGdxia.clear();
635 it->dGdxia.resize(ts2, 0.0);
636 MPI_Unpack(buf, bs, &p, &(it->dGdxia.front()), ts2, MPI_DOUBLE, comm);
637 }
638#endif
639
640 // Atom.dGdr
641 ts2 = 0;
642 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
643 if (ts2 > 0)
644 {
645 it->dGdr.clear();
646 it->dGdr.resize(ts2);
647 for (vector<Vec3D>::iterator it2 = it->dGdr.begin();
648 it2 != it->dGdr.end(); ++it2)
649 {
650 MPI_Unpack(buf, bs, &p, it2->r, 3, MPI_DOUBLE, comm);
651 }
652 }
653
654 // Atom.neighbors
655 ts2 = 0;
656 MPI_Unpack(buf, bs, &p, &ts2, 1, MPI_SIZE_T, comm);
657 if (ts2 > 0)
658 {
659 it->neighbors.clear();
660 it->neighbors.resize(ts2);
661 for (vector<Atom::Neighbor>::iterator it2 =
662 it->neighbors.begin(); it2 != it->neighbors.end(); ++it2)
663 {
664 // Neighbor
665 MPI_Unpack(buf, bs, &p, &(it2->index ), 1, MPI_SIZE_T , comm);
666 MPI_Unpack(buf, bs, &p, &(it2->tag ), 1, MPI_INT64_T, comm);
667 MPI_Unpack(buf, bs, &p, &(it2->element ), 1, MPI_SIZE_T , comm);
668 MPI_Unpack(buf, bs, &p, &(it2->d ), 1, MPI_DOUBLE , comm);
669 MPI_Unpack(buf, bs, &p, it2->dr.r , 3, MPI_DOUBLE , comm);
670
671 size_t ts3 = 0;
672#ifndef N2P2_NO_SF_CACHE
673 // Neighbor.cache
674 ts3 = 0;
675 MPI_Unpack(buf, bs, &p, &ts3, 1, MPI_SIZE_T, comm);
676 if (ts3 > 0)
677 {
678 it2->cache.clear();
679 it2->cache.resize(ts3, 0.0);
680 MPI_Unpack(buf, bs, &p, &(it2->cache.front()), ts3, MPI_DOUBLE, comm);
681 }
682#endif
683
684 // Neighbor.dGdr
685 ts3 = 0;
686 MPI_Unpack(buf, bs, &p, &ts3, 1, MPI_SIZE_T, comm);
687 if (ts3 > 0)
688 {
689 it2->dGdr.clear();
690 it2->dGdr.resize(ts3);
691 for (vector<Vec3D>::iterator it3 = it2->dGdr.begin();
692 it3 != it2->dGdr.end(); ++it3)
693 {
694 MPI_Unpack(buf, bs, &p, it3->r, 3, MPI_DOUBLE, comm);
695 }
696 }
697 }
698 }
699 }
700 }
701
702 delete[] buf;
703
704 return bs;
705}

References nnp::Structure::atoms, nnp::Structure::box, nnp::Structure::charge, nnp::Structure::chargeRef, comm, nnp::Structure::comment, nnp::Structure::energy, nnp::Structure::energyRef, nnp::Structure::hasNeighborList, nnp::Structure::hasSymmetryFunctionDerivatives, nnp::Structure::hasSymmetryFunctions, nnp::Structure::index, nnp::Structure::invbox, nnp::Structure::isPeriodic, nnp::Structure::isTriclinic, MPI_SIZE_T, nnp::Structure::numAtoms, nnp::Structure::numAtomsPerElement, nnp::Structure::numElements, nnp::Structure::numElementsPresent, p, nnp::Structure::pbc, nnp::Vec3D::r, nnp::Structure::sampleType, and nnp::Structure::volume.

Referenced by distributeStructures(), and prepareNumericForces().

Here is the caller graph for this function:

◆ distributeStructures()

int Dataset::distributeStructures ( bool  randomize,
bool  excludeRank0 = false,
std::string const &  fileName = "input.data" 
)

Read data file and distribute structures among processors.

Parameters
[in]randomizeIf true randomly distribute structures, otherwise they are distributed in order.
[in]excludeRank0If true no structures are distributed to MPI process with rank 0.
[in]fileNameData file name, e.g. "input.data".
Returns
Number of bytes transferred.

Definition at line 722 of file Dataset.cpp.

725{
726 log << "\n";
727 log << "*** STRUCTURE DISTRIBUTION **************"
728 "**************************************\n";
729 log << "\n";
730
731 ifstream dataFile;
732 vector<size_t> numStructuresPerProc;
733
734 if (excludeRank0)
735 {
736 log << "No structures will be distributed to rank 0.\n";
737 if (numProcs == 1)
738 {
739 throw runtime_error("ERROR: Can not distribute structures, "
740 "at least 2 MPI tasks are required.\n");
741 }
742 }
743 size_t numReceivers = numProcs;
744 if (excludeRank0) numReceivers--;
745
746 if (myRank == 0)
747 {
748 log << strpr("Reading configurations from data file: %s.\n",
749 fileName.c_str());
750 dataFile.open(fileName.c_str());
752 log << strpr("Total number of structures: %zu\n", numStructures);
753 dataFile.clear();
754 dataFile.seekg(0);
755 }
756 MPI_Bcast(&numStructures, 1, MPI_SIZE_T, 0, comm);
757 if (numStructures < numReceivers)
758 {
759 throw runtime_error("ERROR: More receiving processors than "
760 "structures.\n");
761 }
762
763 numStructuresPerProc.resize(numProcs, 0);
764 if (myRank == 0)
765 {
766 size_t quotient = numStructures / numReceivers;
767 size_t remainder = numStructures % numReceivers;
768 for (size_t i = 0; i < (size_t) numProcs; i++)
769 {
770 if (i != 0 || (!excludeRank0))
771 {
772 numStructuresPerProc.at(i) = quotient;
773 if (remainder > 0 && i > 0 && i <= remainder)
774 {
775 numStructuresPerProc.at(i)++;
776 }
777 }
778 }
779 if (remainder == 0)
780 {
781 log << strpr("Number of structures per processor: %d\n", quotient);
782 }
783 else
784 {
785 log << strpr("Number of structures per"
786 " processor: %d (%d) or %d (%d)\n",
787 quotient,
788 numReceivers - remainder,
789 quotient + 1,
790 remainder);
791 }
792 }
793 MPI_Bcast(&(numStructuresPerProc.front()), numProcs, MPI_SIZE_T, 0, comm);
794
795 int bytesTransferred = 0;
796 size_t numMyStructures = numStructuresPerProc.at(myRank);
797 if (myRank == 0)
798 {
799 size_t countStructures = 0;
800 vector<size_t> countStructuresPerProc;
801
802 countStructuresPerProc.resize(numProcs, 0);
803
804 if (randomize)
805 {
806 while (countStructures < numStructures)
807 {
808 int proc = gsl_rng_uniform_int(rng, numProcs);
809 if (countStructuresPerProc.at(proc)
810 < numStructuresPerProc.at(proc))
811 {
812 if (proc == 0)
813 {
814 structures.push_back(Structure());
815 structures.back().setElementMap(elementMap);
816 structures.back().index = countStructures;
817 structures.back().readFromFile(dataFile);
819 }
820 else
821 {
822 Structure tmpStructure;
823 tmpStructure.setElementMap(elementMap);
824 tmpStructure.index = countStructures;
825 tmpStructure.readFromFile(dataFile);
826 removeEnergyOffset(tmpStructure);
827 bytesTransferred += sendStructure(tmpStructure, proc);
828 }
829 countStructuresPerProc.at(proc)++;
830 countStructures++;
831 }
832 }
833 }
834 else
835 {
836 for (int proc = 0; proc < numProcs; ++proc)
837 {
838 for (size_t i = 0; i < numStructuresPerProc.at(proc); ++i)
839 {
840 if (proc == 0)
841 {
842 structures.push_back(Structure());
843 structures.back().setElementMap(elementMap);
844 structures.back().index = countStructures;
845 structures.back().readFromFile(dataFile);
847 }
848 else
849 {
850 Structure tmpStructure;
851 tmpStructure.setElementMap(elementMap);
852 tmpStructure.index = countStructures;
853 tmpStructure.readFromFile(dataFile);
854 removeEnergyOffset(tmpStructure);
855 bytesTransferred += sendStructure(tmpStructure, proc);
856 }
857 countStructuresPerProc.at(proc)++;
858 countStructures++;
859 }
860 }
861 }
862 dataFile.close();
863 }
864 else
865 {
866 for (size_t i = 0; i < numMyStructures; i++)
867 {
868 structures.push_back(Structure());
869 structures.back().setElementMap(elementMap);
870 bytesTransferred += recvStructure(&(structures.back()), 0);
871 }
872 }
873
874 log << strpr("Distributed %zu structures,"
875 " %d bytes (%.2f MiB) transferred.\n",
877 bytesTransferred,
878 bytesTransferred / 1024. / 1024.);
879 log << strpr("Number of local structures: %zu\n", structures.size());
880 log << "*****************************************"
881 "**************************************\n";
882
883 return bytesTransferred;
884}
int recvStructure(Structure *structure, int src)
Receive one structure from source process.
Definition: Dataset.cpp:462
std::vector< Structure > structures
All structures in this dataset.
Definition: Dataset.h:195
int sendStructure(Structure const &structure, int dest) const
Send one structure to destination process.
Definition: Dataset.cpp:251
std::size_t getNumStructures(std::ifstream &dataFile)
Get number of structures in data file.
Definition: Dataset.cpp:707
ElementMap elementMap
Global element map, populated by setupElementMap().
Definition: Mode.h:508
void removeEnergyOffset(Structure &structure, bool ref=true)
Remove atomic energy offsets from reference energy.
Definition: Mode.cpp:1480
void setElementMap(ElementMap const &elementMap)
Set element map of structure.
Definition: Structure.cpp:60
void readFromFile(std::string const fileName="input.data")
Read configuration from file.
Definition: Structure.cpp:86

References comm, nnp::Mode::elementMap, getNumStructures(), nnp::Structure::index, nnp::Mode::log, MPI_SIZE_T, myRank, numProcs, numStructures, nnp::Structure::readFromFile(), recvStructure(), nnp::Mode::removeEnergyOffset(), rng, sendStructure(), nnp::Structure::setElementMap(), nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prepareNumericForces()

size_t Dataset::prepareNumericForces ( Structure original,
double  delta 
)

Prepare numeric force check for a single structure.

Parameters
[in,out]originalThe structure under investigation.
[in]deltaCentral difference delta for positions
Returns
Number of structures distributed (= 1 + 6 * numAtoms).

Distributes copies of original structure among all processors and modifies positions for central difference.

Definition at line 886 of file Dataset.cpp.

887{
888 // Be sure that structure memory is cleared.
889 vector<Structure>().swap(structures);
890
891 // Send original to all processors.
892 for (int p = 1; p < numProcs; ++p)
893 {
894 if (myRank == 0) sendStructure(original, p);
895 else if (myRank == p) recvStructure(&original, 0);
896 }
897
898 vector<size_t> numStructuresPerProc(numProcs, 0);
899 // Number of structures that need to be calculated:
900 // - Original structure +
901 // - Two times for each force component.
902 numStructures = 1 + 6 * original.numAtoms;
903 size_t quotient = numStructures / numProcs;
904 size_t remainder = numStructures % numProcs;
905 for (size_t i = 0; i < (size_t) numProcs; i++)
906 {
907 numStructuresPerProc.at(i) = quotient;
908 if (remainder > 0 && i < remainder) numStructuresPerProc.at(i)++;
909 }
910
911 // Rank 0 needs an unaltered version.
912 if (myRank == 0)
913 {
914 structures.push_back(original);
915 structures.back().setElementMap(elementMap);
916 structures.back().comment = "original";
917 }
918 // Now make as many copies as each processor needs.
919 size_t count = 0;
920 size_t iAtom = 0;
921 size_t ixyz = 0;
922 int sign = 0;
923 for (int p = 0; p < numProcs; ++p)
924 {
925 size_t istart = 0;
926 if (p == 0) istart = 1;
927 for (size_t i = istart; i < numStructuresPerProc.at(p); ++i)
928 {
929 iAtom = count / 6;
930 ixyz = count % 3;
931 sign = 2 * ((count % 6) / 3) - 1;
932 if (p == myRank)
933 {
934 structures.push_back(original);
935 Structure& s = structures.back();
937 // Write identifiers to comment field.
938 s.comment = strpr("%zu %zu %zu %d", count, iAtom, ixyz, sign);
939 // Modify atom position.
940 s.atoms.at(iAtom).r[ixyz] += sign * delta;
941 if (s.isPeriodic) s.remap(s.atoms.at(iAtom));
942 }
943 count++;
944 }
945 }
946
947 return numStructures;
948}
void remap()
Translate all atoms back into box if outside.
Definition: Structure.cpp:443

References nnp::Structure::atoms, nnp::Structure::comment, nnp::Mode::elementMap, nnp::Structure::isPeriodic, myRank, nnp::Structure::numAtoms, numProcs, numStructures, p, recvStructure(), nnp::Structure::remap(), sendStructure(), nnp::Structure::setElementMap(), nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ toNormalizedUnits()

void Dataset::toNormalizedUnits ( )

Switch all structures to normalized units.

Definition at line 950 of file Dataset.cpp.

951{
952 for (vector<Structure>::iterator it = structures.begin();
953 it != structures.end(); ++it)
954 {
955 it->toNormalizedUnits(meanEnergy, convEnergy, convLength);
956 }
957
958 return;
959}
double convEnergy
Definition: Mode.h:523
double convLength
Definition: Mode.h:524
double meanEnergy
Definition: Mode.h:522

References nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Mode::meanEnergy, and structures.

Referenced by main().

Here is the caller graph for this function:

◆ toPhysicalUnits()

void Dataset::toPhysicalUnits ( )

Switch all structures to physical units.

Definition at line 961 of file Dataset.cpp.

962{
963 for (vector<Structure>::iterator it = structures.begin();
964 it != structures.end(); ++it)
965 {
966 it->toPhysicalUnits(meanEnergy, convEnergy, convLength);
967 }
968
969 return;
970}

References nnp::Mode::convEnergy, nnp::Mode::convLength, nnp::Mode::meanEnergy, and structures.

Referenced by main().

Here is the caller graph for this function:

◆ collectSymmetryFunctionStatistics()

void Dataset::collectSymmetryFunctionStatistics ( )

Collect symmetry function statistics from all processors.

Definition at line 972 of file Dataset.cpp.

973{
974 for (vector<Element>::iterator it = elements.begin();
975 it != elements.end(); ++it)
976 {
977 // If no atoms of this element exist on this proc, create empty
978 // statistics.
979 if (it->statistics.data.size() == 0)
980 {
981 log << strpr("WARNING: No statistics for element %zu (%2s) found, "
982 "process %d has no corresponding atoms, creating "
983 "empty statistics.\n",
984 it->getIndex(),
985 it->getSymbol().c_str(),
986 myRank);
987 }
988 for (size_t i = 0; i < it->numSymmetryFunctions(); ++i)
989 {
990 SymFncStatistics::Container& c = it->statistics.data[i];
991 MPI_Allreduce(MPI_IN_PLACE, &(c.count), 1, MPI_SIZE_T, MPI_SUM, comm);
992 MPI_Allreduce(MPI_IN_PLACE, &(c.min ), 1, MPI_DOUBLE, MPI_MIN, comm);
993 MPI_Allreduce(MPI_IN_PLACE, &(c.max ), 1, MPI_DOUBLE, MPI_MAX, comm);
994 MPI_Allreduce(MPI_IN_PLACE, &(c.sum ), 1, MPI_DOUBLE, MPI_SUM, comm);
995 MPI_Allreduce(MPI_IN_PLACE, &(c.sum2 ), 1, MPI_DOUBLE, MPI_SUM, comm);
996 }
997 }
998
999 return;
1000}
std::vector< Element > elements
Definition: Mode.h:528
Struct containing statistics gathered during symmetry function calculation.
double min
Minimum symmetry function value encountered.
double sum2
Sum of squared symmetry function values (to compute sigma).
double sum
Sum of symmetry function values (to compute mean).
double max
Maximum symmetry function value encountered.
std::size_t count
Counts total number of symmetry function evaluations.

References comm, nnp::SymFncStatistics::Container::count, nnp::Mode::elements, nnp::Mode::log, nnp::SymFncStatistics::Container::max, nnp::SymFncStatistics::Container::min, MPI_SIZE_T, myRank, nnp::strpr(), nnp::SymFncStatistics::Container::sum, and nnp::SymFncStatistics::Container::sum2.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeSymmetryFunctionScaling()

void Dataset::writeSymmetryFunctionScaling ( std::string const &  fileName = "scaling.data")

Write symmetry function scaling values to file.

Parameters
[in]fileNameScaling data file name (e.g. "scaling.data").

Definition at line 1002 of file Dataset.cpp.

1003{
1004 log << "\n";
1005 log << "*** SYMMETRY FUNCTION SCALING ***********"
1006 "**************************************\n";
1007 log << "\n";
1008
1009 if (myRank == 0)
1010 {
1011 log << strpr("Writing symmetry function scaling file: %s.\n",
1012 fileName.c_str());
1013 ofstream sFile;
1014 sFile.open(fileName.c_str());
1015
1016 // File header.
1017 vector<string> title;
1018 vector<string> colName;
1019 vector<string> colInfo;
1020 vector<size_t> colSize;
1021 title.push_back("Symmetry function scaling data.");
1022 colSize.push_back(10);
1023 colName.push_back("e_index");
1024 colInfo.push_back("Element index.");
1025 colSize.push_back(10);
1026 colName.push_back("sf_index");
1027 colInfo.push_back("Symmetry function index.");
1028 colSize.push_back(24);
1029 colName.push_back("sf_min");
1030 colInfo.push_back("Symmetry function minimum.");
1031 colSize.push_back(24);
1032 colName.push_back("sf_max");
1033 colInfo.push_back("Symmetry function maximum.");
1034 colSize.push_back(24);
1035 colName.push_back("sf_mean");
1036 colInfo.push_back("Symmetry function mean.");
1037 colSize.push_back(24);
1038 colName.push_back("sf_sigma");
1039 colInfo.push_back("Symmetry function sigma.");
1040 appendLinesToFile(sFile,
1041 createFileHeader(title, colSize, colName, colInfo));
1042
1043 log << "\n";
1044 log << "Abbreviations:\n";
1045 log << "--------------\n";
1046 log << "ind ...... Symmetry function index.\n";
1047 log << "min ...... Minimum symmetry function value.\n";
1048 log << "max ...... Maximum symmetry function value.\n";
1049 log << "mean ..... Mean symmetry function value.\n";
1050 log << "sigma .... Standard deviation of symmetry function values.\n";
1051 log << "spread ... (max - min) / sigma.\n";
1052 log << "\n";
1053 for (vector<Element>::const_iterator it = elements.begin();
1054 it != elements.end(); ++it)
1055 {
1056 log << strpr("Scaling data for symmetry functions element %2s :\n",
1057 it->getSymbol().c_str());
1058 log << "-----------------------------------------"
1059 "--------------------------------------\n";
1060 log << " ind min max mean sigma spread\n";
1061 log << "-----------------------------------------"
1062 "--------------------------------------\n";
1063 for (size_t i = 0; i < it->numSymmetryFunctions(); ++i)
1064 {
1066 = it->statistics.data.at(i);
1067 size_t n = c.count;
1068 double const mean = c.sum / n;
1069 double const sigma = sqrt((c.sum2 - c.sum * c.sum / n)
1070 / (n - 1));
1071 double const spread = (c.max - c.min) / sigma;
1072 sFile << strpr("%10d %10d %24.16E %24.16E %24.16E %24.16E\n",
1073 it->getIndex() + 1,
1074 i + 1,
1075 c.min,
1076 c.max,
1077 mean,
1078 sigma);
1079 log << strpr("%4zu %9.2E %9.2E %9.2E %9.2E %9.2E\n",
1080 i + 1,
1081 c.min,
1082 c.max,
1083 mean,
1084 sigma,
1085 spread);
1086 }
1087 log << "-----------------------------------------"
1088 "--------------------------------------\n";
1089 }
1090 // Finally decided to remove this legacy line...
1091 //sFile << strpr("%f %f\n", 0.0, 0.0);
1092 sFile.close();
1093 }
1094
1095 log << "*****************************************"
1096 "**************************************\n";
1097
1098 return;
1099}
vector< string > createFileHeader(vector< string > const &title, vector< size_t > const &colSize, vector< string > const &colName, vector< string > const &colInfo, char const &commentChar)
Definition: utility.cpp:104
void appendLinesToFile(ofstream &file, vector< string > const lines)
Append multiple lines of strings to open file stream.
Definition: utility.cpp:218

References nnp::appendLinesToFile(), nnp::SymFncStatistics::Container::count, nnp::createFileHeader(), nnp::Mode::elements, nnp::Mode::log, nnp::SymFncStatistics::Container::max, nnp::SymFncStatistics::Container::min, myRank, nnp::strpr(), nnp::SymFncStatistics::Container::sum, and nnp::SymFncStatistics::Container::sum2.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeSymmetryFunctionHistograms()

void Dataset::writeSymmetryFunctionHistograms ( std::size_t  numBins,
std::string  fileNameFormat = "sf.%03zu.%04zu.histo" 
)

Calculate and write symmetry function histograms.

Parameters
[in]numBinsNumber of histogram bins used.
[in]fileNameFormatFormat for histogram file names, must include placeholders for element and symmetry function number.

Definition at line 1101 of file Dataset.cpp.

1103{
1104 log << "\n";
1105 log << "*** SYMMETRY FUNCTION HISTOGRAMS ********"
1106 "**************************************\n";
1107 log << "\n";
1108
1109 // Initialize histograms.
1110 numBins--;
1111 vector<vector<gsl_histogram*> > histograms;
1112 for (vector<Element>::const_iterator it = elements.begin();
1113 it != elements.end(); ++it)
1114 {
1115 histograms.push_back(vector<gsl_histogram*>());
1116 for (size_t i = 0; i < it->numSymmetryFunctions(); ++i)
1117 {
1118 double l = safeFind(it->statistics.data, i).min;
1119 double h = safeFind(it->statistics.data, i).max;
1120 if (l < h)
1121 {
1122 // Add an extra bin at the end to cover complete range.
1123 h += (h - l) / numBins;
1124 histograms.back().push_back(gsl_histogram_alloc(numBins + 1));
1125 gsl_histogram_set_ranges_uniform(histograms.back().back(),
1126 l,
1127 h);
1128 }
1129 else
1130 {
1131 // Use nullptr so signalize non-existing histogram.
1132 histograms.back().push_back(nullptr);
1133 log << strpr("WARNING: Symmetry function min equals max, "
1134 "ommitting histogram (Element %2s SF %4zu "
1135 "(line %4zu).\n",
1136 it->getSymbol().c_str(),
1137 i,
1138 it->getSymmetryFunction(i).getLineNumber() + 1);
1139 }
1140 }
1141 }
1142
1143 // Increment histograms with symmetry function values.
1144 for (vector<Structure>::const_iterator it = structures.begin();
1145 it != structures.end(); ++it)
1146 {
1147 for (vector<Atom>::const_iterator it2 = it->atoms.begin();
1148 it2 != it->atoms.end(); ++it2)
1149 {
1150 size_t e = it2->element;
1151 vector<gsl_histogram*>& h = histograms.at(e);
1152 for (size_t s = 0; s < it2->G.size(); ++s)
1153 {
1154 if (h.at(s) == nullptr) continue;
1155 gsl_histogram_increment(h.at(s), it2->G.at(s));
1156 }
1157 }
1158 }
1159
1160 // Collect histograms from all processes.
1161 for (vector<vector<gsl_histogram*> >::const_iterator it =
1162 histograms.begin(); it != histograms.end(); ++it)
1163 {
1164 for (vector<gsl_histogram*>::const_iterator it2 = it->begin();
1165 it2 != it->end(); ++it2)
1166 {
1167 if ((*it2) == nullptr) continue;
1168 MPI_Allreduce(MPI_IN_PLACE, (*it2)->bin, numBins + 1, MPI_DOUBLE, MPI_SUM, comm);
1169 }
1170 }
1171
1172 // Write histogram to file.
1173 if (myRank == 0)
1174 {
1175 log << strpr("Writing histograms with %zu bins.\n", numBins + 1);
1176 for (size_t e = 0; e < elements.size(); ++e)
1177 {
1178 for (size_t s = 0; s < elements.at(e).numSymmetryFunctions(); ++s)
1179 {
1180 gsl_histogram*& h = histograms.at(e).at(s);
1181 if (h == nullptr) continue;
1182 FILE* fp = 0;
1183 string fileName = strpr(fileNameFormat.c_str(),
1185 s + 1);
1186 fp = fopen(fileName.c_str(), "w");
1187 if (fp == 0)
1188 {
1189 throw runtime_error(strpr("ERROR: Could not open file:"
1190 " %s.\n", fileName.c_str()));
1191 }
1192 vector<string> info = elements.at(e).infoSymmetryFunction(s);
1193 for (vector<string>::const_iterator it = info.begin();
1194 it != info.end(); ++it)
1195 {
1196 fprintf(fp, "#SFINFO %s\n", it->c_str());
1197 }
1199 = elements.at(e).statistics.data.at(s);
1200 size_t n = c.count;
1201 fprintf(fp, "#SFINFO min %15.8E\n", c.min);
1202 fprintf(fp, "#SFINFO max %15.8E\n", c.max);
1203 fprintf(fp, "#SFINFO mean %15.8E\n", c.sum / n);
1204 fprintf(fp, "#SFINFO sigma %15.8E\n",
1205 sqrt(1.0 / (n - 1) * (c.sum2 - c.sum * c.sum / n)));
1206
1207 // File header.
1208 vector<string> title;
1209 vector<string> colName;
1210 vector<string> colInfo;
1211 vector<size_t> colSize;
1212 title.push_back("Symmetry function histogram.");
1213 colSize.push_back(16);
1214 colName.push_back("symfunc_l");
1215 colInfo.push_back("Symmetry function value, left bin limit.");
1216 colSize.push_back(16);
1217 colName.push_back("symfunc_r");
1218 colInfo.push_back("Symmetry function value, right bin limit.");
1219 colSize.push_back(16);
1220 colName.push_back("hist");
1221 colInfo.push_back("Histogram count.");
1223 createFileHeader(title,
1224 colSize,
1225 colName,
1226 colInfo));
1227
1228 gsl_histogram_fprintf(fp, h, "%16.8E", "%16.8E");
1229 fflush(fp);
1230 fclose(fp);
1231 fp = 0;
1232 }
1233 }
1234 }
1235
1236 for (vector<vector<gsl_histogram*> >::const_iterator it =
1237 histograms.begin(); it != histograms.end(); ++it)
1238 {
1239 for (vector<gsl_histogram*>::const_iterator it2 = it->begin();
1240 it2 != it->end(); ++it2)
1241 {
1242 if ((*it2) == nullptr) continue;
1243 gsl_histogram_free(*it2);
1244 }
1245 }
1246
1247 log << "*****************************************"
1248 "**************************************\n";
1249
1250 return;
1251}
std::size_t atomicNumber(std::size_t index) const
Get atomic number from element index.
Definition: ElementMap.h:145
V const & safeFind(std::map< K, V > const &stdMap, typename std::map< K, V >::key_type const &key)
Safely access map entry.
Definition: utility.h:135

References nnp::appendLinesToFile(), nnp::ElementMap::atomicNumber(), comm, nnp::SymFncStatistics::Container::count, nnp::createFileHeader(), nnp::Mode::elementMap, nnp::Mode::elements, nnp::Mode::log, nnp::SymFncStatistics::Container::max, nnp::SymFncStatistics::Container::min, myRank, nnp::safeFind(), nnp::strpr(), structures, nnp::SymFncStatistics::Container::sum, and nnp::SymFncStatistics::Container::sum2.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeSymmetryFunctionFile()

void Dataset::writeSymmetryFunctionFile ( std::string  fileName = "function.data")

Write symmetry function legacy file ("function.data").

Parameters
[in]fileNameFile name for symmetry function file.

Definition at line 1253 of file Dataset.cpp.

1254{
1255 log << "\n";
1256 log << "*** SYMMETRY FUNCTION FILE **************"
1257 "**************************************\n";
1258 log << "\n";
1259
1260 // Create empty file.
1261 log << strpr("Writing symmetry functions to file: %s\n", fileName.c_str());
1262 if (myRank == 0)
1263 {
1264 ofstream file;
1265 file.open(fileName.c_str());
1266 file.close();
1267 }
1268 MPI_Barrier(comm);
1269
1270 // Prepare structure iterator.
1271 vector<Structure>::const_iterator it = structures.begin();
1272 // Loop over all structures (on each proc the local structures are stored
1273 // with increasing index).
1274 for (size_t i = 0; i < numStructures; ++i)
1275 {
1276 // If this proc holds structure with matching index,
1277 // open file and write symmetry functions.
1278 if (it != structures.end() && i == it->index)
1279 {
1280 ofstream file;
1281 file.open(fileName.c_str(), ios_base::app);
1282 file << strpr("%6zu\n", it->numAtoms);
1283 // Loop over atoms.
1284 for (vector<Atom>::const_iterator it2 = it->atoms.begin();
1285 it2 != it->atoms.end(); ++it2)
1286 {
1287 // Loop over symmetry functions.
1288 file << strpr("%3zu ", elementMap.atomicNumber(it2->element));
1289 for (vector<double>::const_iterator it3 = it2->G.begin();
1290 it3 != it2->G.end(); ++it3)
1291 {
1292 file << strpr(" %14.10f", *it3);
1293 }
1294 file << '\n';
1295 }
1296 // There is no charge NN, so first and last entry is zero.
1297 double energy = 0.0;
1298 if (normalize) energy = physicalEnergy(*it);
1299 else energy = it->energyRef;
1300 energy += getEnergyOffset(*it);
1301 energy /= it->numAtoms;
1302 file << strpr(" %19.10f %19.10f %19.10f %19.10f\n",
1303 0.0, energy, energy, 0.0);
1304 file.flush();
1305 file.close();
1306 // Iterate to next structure.
1307 ++it;
1308 }
1309 MPI_Barrier(comm);
1310 }
1311
1312 log << "*****************************************"
1313 "**************************************\n";
1314
1315 return;
1316}
double physicalEnergy(Structure const &structure, bool ref=true) const
Undo normalization for a given energy of structure.
Definition: Mode.cpp:1557
bool normalize
Definition: Mode.h:514
double getEnergyOffset(Structure const &structure) const
Get atomic energy offset for given structure.
Definition: Mode.cpp:1499

References nnp::ElementMap::atomicNumber(), comm, nnp::Mode::elementMap, nnp::Mode::getEnergyOffset(), nnp::Mode::log, myRank, nnp::Mode::normalize, numStructures, nnp::Mode::physicalEnergy(), nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeNeighborHistogram()

size_t Dataset::writeNeighborHistogram ( std::string const &  fileNameHisto = "neighbors.histo",
std::string const &  fileNameStructure = "neighbors.out" 
)

Calculate and write neighbor histogram and per-structure statistics.

Parameters
[in]fileNameHistoName of histogram file.
[in]fileNameStructureName of per-structure file.
Returns
Maximum number of neighbors.

Definition at line 1318 of file Dataset.cpp.

1320{
1321 log << "\n";
1322 log << "*** NEIGHBOR STATISTICS/HISTOGRAM *******"
1323 "**************************************\n";
1324 log << "\n";
1325
1326 ofstream statFile;
1327 string fileNameLocal = strpr("%s.%04d", fileNameStructure.c_str(), myRank);
1328 statFile.open(fileNameLocal.c_str());
1329 if (myRank == 0)
1330 {
1331 // File header.
1332 vector<string> title;
1333 vector<string> colName;
1334 vector<string> colInfo;
1335 vector<size_t> colSize;
1336 title.push_back("Neighbor statistics for each structure.");
1337 colSize.push_back(10);
1338 colName.push_back("struct");
1339 colInfo.push_back("Index of structure (starting with 1).");
1340 colSize.push_back(10);
1341 colName.push_back("natoms");
1342 colInfo.push_back("Number of atoms in structure.");
1343 colSize.push_back(10);
1344 colName.push_back("min");
1345 colInfo.push_back("Minimum number of neighbors over all atoms.");
1346 colSize.push_back(10);
1347 colName.push_back("max");
1348 colInfo.push_back("Maximum number of neighbors over all atoms.");
1349 colSize.push_back(16);
1350 colName.push_back("mean");
1351 colInfo.push_back("Mean number of neighbors over all atoms.");
1352 appendLinesToFile(statFile,
1353 createFileHeader(title, colSize, colName, colInfo));
1354 }
1355
1356 // Compute per-structure statistics and prepare histogram boundaries.
1357 size_t numAtomsGlobal = 0;
1358 size_t minNeighborsGlobal = numeric_limits<size_t>::max();
1359 size_t maxNeighborsGlobal = 0;
1360 double meanNeighborsGlobal = 0.0;
1361 for (auto const& s : structures)
1362 {
1363 size_t minNeighbors = numeric_limits<size_t>::max();
1364 size_t maxNeighbors = 0;
1365 double meanNeighbors = 0.0;
1366 for (auto const& a : s.atoms)
1367 {
1368 size_t const n = a.numNeighbors;
1369 minNeighbors = min(minNeighbors, n);
1370 maxNeighbors = max(maxNeighbors, n);
1371 meanNeighbors += n;
1372 }
1373 numAtomsGlobal += s.numAtoms;
1374 minNeighborsGlobal = min(minNeighborsGlobal, minNeighbors);
1375 maxNeighborsGlobal = max(maxNeighborsGlobal, maxNeighbors);
1376 meanNeighborsGlobal += meanNeighbors;
1377 statFile << strpr("%10zu %10zu %10zu %10zu %16.8E\n",
1378 s.index + 1,
1379 s.numAtoms,
1381 maxNeighbors,
1382 meanNeighbors / s.numAtoms);
1383 }
1384 statFile.flush();
1385 statFile.close();
1386 MPI_Allreduce(MPI_IN_PLACE, &numAtomsGlobal , 1, MPI_SIZE_T, MPI_SUM, comm);
1387 MPI_Allreduce(MPI_IN_PLACE, &minNeighborsGlobal , 1, MPI_SIZE_T, MPI_MIN, comm);
1388 MPI_Allreduce(MPI_IN_PLACE, &maxNeighborsGlobal , 1, MPI_SIZE_T, MPI_MAX, comm);
1389 MPI_Allreduce(MPI_IN_PLACE, &meanNeighborsGlobal, 1, MPI_DOUBLE, MPI_SUM, comm);
1390 log << strpr("Minimum number of neighbors: %d\n", minNeighborsGlobal);
1391 log << strpr("Mean number of neighbors: %.1f\n",
1392 meanNeighborsGlobal / numAtomsGlobal);
1393 log << strpr("Maximum number of neighbors: %d\n", maxNeighborsGlobal);
1394
1395 // Calculate neighbor histogram.
1396 gsl_histogram* histNeighbors = NULL;
1397 histNeighbors = gsl_histogram_alloc(maxNeighborsGlobal + 1);
1398 gsl_histogram_set_ranges_uniform(histNeighbors,
1399 -0.5,
1400 maxNeighborsGlobal + 0.5);
1401 for (vector<Structure>::const_iterator it = structures.begin();
1402 it != structures.end(); ++it)
1403 {
1404 for (vector<Atom>::const_iterator it2 = it->atoms.begin();
1405 it2 != it->atoms.end(); ++it2)
1406 {
1407 gsl_histogram_increment(histNeighbors, it2->numNeighbors);
1408 }
1409 }
1410 MPI_Allreduce(MPI_IN_PLACE, histNeighbors->bin, maxNeighborsGlobal + 1, MPI_DOUBLE, MPI_SUM, comm);
1411
1412 // Write histogram to file.
1413 if (myRank == 0)
1414 {
1415 log << strpr("Neighbor histogram file: %s.\n", fileNameHisto.c_str());
1416 FILE* fp = 0;
1417 fp = fopen(fileNameHisto.c_str(), "w");
1418 if (fp == 0)
1419 {
1420 throw runtime_error(strpr("ERROR: Could not open file: %s.\n",
1421 fileNameHisto.c_str()));
1422 }
1423
1424 // File header.
1425 vector<string> title;
1426 vector<string> colName;
1427 vector<string> colInfo;
1428 vector<size_t> colSize;
1429 title.push_back("Neighbor count histogram.");
1430 colSize.push_back(9);
1431 colName.push_back("neigh_l");
1432 colInfo.push_back("Number of neighbors, left bin limit.");
1433 colSize.push_back(9);
1434 colName.push_back("neigh_r");
1435 colInfo.push_back("Number of neighbors, right bin limit.");
1436 colSize.push_back(16);
1437 colName.push_back("hist");
1438 colInfo.push_back("Histogram count.");
1440 createFileHeader(title, colSize, colName, colInfo));
1441
1442 gsl_histogram_fprintf(fp, histNeighbors, "%9.1f", "%16.8E");
1443 fflush(fp);
1444 fclose(fp);
1445 fp = 0;
1446 }
1447
1448 MPI_Barrier(comm);
1449 if (myRank == 0)
1450 {
1451 log << strpr("Combining per-structure neighbor statistics file: %s.\n",
1452 fileNameStructure.c_str());
1453 combineFiles(fileNameStructure);
1454 }
1455
1456 gsl_histogram_free(histNeighbors);
1457
1458 log << "*****************************************"
1459 "**************************************\n";
1460
1461 return maxNeighborsGlobal;
1462}
void combineFiles(std::string filePrefix) const
Combine individual MPI proc files to one.
Definition: Dataset.cpp:1742
std::vector< std::size_t > minNeighbors
Definition: Mode.h:518

References nnp::appendLinesToFile(), combineFiles(), comm, nnp::createFileHeader(), nnp::Mode::log, nnp::Mode::minNeighbors, MPI_SIZE_T, myRank, nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sortNeighborLists()

void Dataset::sortNeighborLists ( )

Sort all neighbor lists according to element and distance.

Definition at line 1464 of file Dataset.cpp.

1465{
1466 log << "\n";
1467 log << "*** NEIGHBOR LIST ***********************"
1468 "**************************************\n";
1469 log << "\n";
1470
1471 log << "Sorting neighbor lists according to element and distance.\n";
1472
1473 for (vector<Structure>::iterator it = structures.begin();
1474 it != structures.end(); ++it)
1475 {
1476 for (vector<Atom>::iterator it2 = it->atoms.begin();
1477 it2 != it->atoms.end(); ++it2)
1478 {
1479 sort(it2->neighbors.begin(), it2->neighbors.end());
1480 }
1481 }
1482
1483 log << "*****************************************"
1484 "**************************************\n";
1485
1486 return;
1487}

References nnp::Mode::log, and structures.

Referenced by main().

Here is the caller graph for this function:

◆ writeNeighborLists()

void Dataset::writeNeighborLists ( std::string const &  fileName = "neighbor-list.data")

Write neighbor list file.

Parameters
[in]fileNameName for neighbor list file.

Definition at line 1488 of file Dataset.cpp.

1489{
1490 log << "\n";
1491 log << "*** NEIGHBOR LIST ***********************"
1492 "**************************************\n";
1493 log << "\n";
1494
1495 string fileNameLocal = strpr("%s.%04d", fileName.c_str(), myRank);
1496 ofstream fileLocal;
1497 fileLocal.open(fileNameLocal.c_str());
1498
1499 for (vector<Structure>::const_iterator it = structures.begin();
1500 it != structures.end(); ++it)
1501 {
1502 fileLocal << strpr("%zu\n", it->numAtoms);
1503 for (vector<Atom>::const_iterator it2 = it->atoms.begin();
1504 it2 != it->atoms.end(); ++it2)
1505 {
1506 fileLocal << strpr("%zu", elementMap.atomicNumber(it2->element));
1507 for (size_t i = 0; i < numElements; ++i)
1508 {
1509 fileLocal << strpr(" %zu", it2->numNeighborsPerElement.at(i));
1510 }
1511 for (vector<Atom::Neighbor>::const_iterator it3
1512 = it2->neighbors.begin(); it3 != it2->neighbors.end(); ++it3)
1513 {
1514 fileLocal << strpr(" %zu", it3->index);
1515 }
1516 fileLocal << '\n';
1517 }
1518 }
1519
1520 fileLocal.flush();
1521 fileLocal.close();
1522 MPI_Barrier(comm);
1523
1524 log << strpr("Writing neighbor lists to file: %s.\n", fileName.c_str());
1525
1526 if (myRank == 0) combineFiles(fileName);
1527
1528 log << "*****************************************"
1529 "**************************************\n";
1530
1531 return;
1532}
std::size_t numElements
Definition: Mode.h:517

References nnp::ElementMap::atomicNumber(), combineFiles(), comm, nnp::Mode::elementMap, nnp::Mode::log, myRank, nnp::Mode::numElements, nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeAtomicEnvironmentFile()

void Dataset::writeAtomicEnvironmentFile ( std::vector< std::vector< std::size_t > >  neighCutoff,
bool  derivatives,
std::string const &  fileNamePrefix = "atomic-env" 
)

Write atomic environment file.

Parameters
[in]neighCutoffMaximum number of neighbor to consider (for each element combination).
[in]derivativesIf true, write separate files for derivates.
[in]fileNamePrefixPrefix for atomic environment files.

This file is used for symmetry function clustering analysis.

Definition at line 1534 of file Dataset.cpp.

1538{
1539 log << "\n";
1540 log << "*** ATOMIC ENVIRONMENT ******************"
1541 "**************************************\n";
1542 log << "\n";
1543
1544 string const fileNamePrefixG = strpr("%s.G" , fileNamePrefix.c_str());
1545 string const fileNamePrefixdGdx = strpr("%s.dGdx", fileNamePrefix.c_str());
1546 string const fileNamePrefixdGdy = strpr("%s.dGdy", fileNamePrefix.c_str());
1547 string const fileNamePrefixdGdz = strpr("%s.dGdz", fileNamePrefix.c_str());
1548
1549 string const fileNameLocalG = strpr("%s.%04d",
1550 fileNamePrefixG.c_str(), myRank);
1551 string const fileNameLocaldGdx = strpr("%s.%04d",
1552 fileNamePrefixdGdx.c_str(), myRank);
1553 string const fileNameLocaldGdy = strpr("%s.%04d",
1554 fileNamePrefixdGdy.c_str(), myRank);
1555 string const fileNameLocaldGdz = strpr("%s.%04d",
1556 fileNamePrefixdGdz.c_str(), myRank);
1557
1558 ofstream fileLocalG;
1559 ofstream fileLocaldGdx;
1560 ofstream fileLocaldGdy;
1561 ofstream fileLocaldGdz;
1562
1563 fileLocalG.open(fileNameLocalG.c_str());
1564 if (derivatives)
1565 {
1566 fileLocaldGdx.open(fileNameLocaldGdx.c_str());
1567 fileLocaldGdy.open(fileNameLocaldGdy.c_str());
1568 fileLocaldGdz.open(fileNameLocaldGdz.c_str());
1569 }
1570
1571 log << "Preparing symmetry functions for atomic environment file(s).\n";
1572 for (size_t i = 0; i < numElements; ++i)
1573 {
1574 for (size_t j = 0; j < numElements; ++j)
1575 {
1576 log << strpr("Maximum number of %2s neighbors for central %2s "
1577 "atoms: %zu\n",
1578 elementMap[j].c_str(),
1579 elementMap[i].c_str(),
1580 neighCutoff.at(i).at(j));
1581 }
1582 }
1583
1584 vector<size_t> neighCount(numElements, 0);
1585 for (vector<Structure>::const_iterator its = structures.begin();
1586 its != structures.end(); ++its)
1587 {
1588 for (vector<Atom>::const_iterator ita = its->atoms.begin();
1589 ita != its->atoms.end(); ++ita)
1590 {
1591 size_t const ea = ita->element;
1592 for (size_t i = 0; i < numElements; ++i)
1593 {
1594 if (ita->numNeighborsPerElement.at(i)
1595 < neighCutoff.at(ea).at(i))
1596 {
1597 throw runtime_error(strpr(
1598 "ERROR: Not enough neighbor atoms, cannot create "
1599 "atomic environment file. Reduce neighbor cutoff for "
1600 "element %2s.\n", elementMap[i].c_str()).c_str());
1601 }
1602 }
1603 fileLocalG << strpr("%2s", elementMap[ita->element].c_str());
1604 fileLocaldGdx << strpr("%2s", elementMap[ita->element].c_str());
1605 fileLocaldGdy << strpr("%2s", elementMap[ita->element].c_str());
1606 fileLocaldGdz << strpr("%2s", elementMap[ita->element].c_str());
1607 // Write atom's own symmetry functions (and derivatives).
1608 for (vector<double>::const_iterator it = ita->G.begin();
1609 it != ita->G.end(); ++it)
1610 {
1611 fileLocalG << strpr(" %16.8E", (*it));
1612 }
1613 if (derivatives)
1614 {
1615 for (vector<Vec3D>::const_iterator it = ita->dGdr.begin();
1616 it != ita->dGdr.end(); ++it)
1617 {
1618 fileLocaldGdx << strpr(" %16.8E", (*it)[0]);
1619 fileLocaldGdy << strpr(" %16.8E", (*it)[1]);
1620 fileLocaldGdz << strpr(" %16.8E", (*it)[2]);
1621 }
1622 }
1623 // Write symmetry functions of neighbors
1624 for (vector<Atom::Neighbor>::const_iterator itn
1625 = ita->neighbors.begin(); itn != ita->neighbors.end(); ++itn)
1626 {
1627 size_t const i = itn->index;
1628 size_t const en = itn->element;
1629 if (neighCount.at(en) < neighCutoff.at(ea).at(en))
1630 {
1631 // Look up symmetry function at Atom instance of neighbor.
1632 Atom const& a = its->atoms.at(i);
1633 for (vector<double>::const_iterator it = a.G.begin();
1634 it != a.G.end(); ++it)
1635 {
1636 fileLocalG << strpr(" %16.8E", (*it));
1637 }
1638 // Log derivatives directly from Neighbor instance.
1639 if (derivatives)
1640 {
1641 // Find atom in neighbor list of neighbor atom.
1642 vector<Atom::Neighbor>::const_iterator itan = find_if(
1643 a.neighbors.begin(), a.neighbors.end(),
1644 [&ita](Atom::Neighbor const& n)
1645 {
1646 return n.index == ita->index;
1647 });
1648 for (vector<Vec3D>::const_iterator it
1649 = itan->dGdr.begin();
1650 it != itan->dGdr.end(); ++it)
1651 {
1652 fileLocaldGdx << strpr(" %16.8E", (*it)[0]);
1653 fileLocaldGdy << strpr(" %16.8E", (*it)[1]);
1654 fileLocaldGdz << strpr(" %16.8E", (*it)[2]);
1655 }
1656 }
1657 neighCount.at(en)++;
1658 }
1659 }
1660 fileLocalG << '\n';
1661 if (derivatives)
1662 {
1663 fileLocaldGdx << '\n';
1664 fileLocaldGdy << '\n';
1665 fileLocaldGdz << '\n';
1666 }
1667 // Reset neighbor counter.
1668 fill(neighCount.begin(), neighCount.end(), 0);
1669 }
1670 }
1671
1672 fileLocalG.flush();
1673 fileLocalG.close();
1674 if (derivatives)
1675 {
1676 fileLocaldGdx.flush();
1677 fileLocaldGdx.close();
1678 fileLocaldGdy.flush();
1679 fileLocaldGdy.close();
1680 fileLocaldGdz.flush();
1681 fileLocaldGdz.close();
1682 }
1683 MPI_Barrier(comm);
1684
1685 if (myRank == 0)
1686 {
1687 log << strpr("Combining atomic environment file: %s.\n",
1688 fileNamePrefixG.c_str());
1689 combineFiles(fileNamePrefixG);
1690 if (derivatives)
1691 {
1692 log << strpr("Combining atomic environment file: %s.\n",
1693 fileNamePrefixdGdx.c_str());
1694 combineFiles(fileNamePrefixdGdx);
1695 log << strpr("Combining atomic environment file: %s.\n",
1696 fileNamePrefixdGdy.c_str());
1697 combineFiles(fileNamePrefixdGdy);
1698 log << strpr("Combining atomic environment file: %s.\n",
1699 fileNamePrefixdGdz.c_str());
1700 combineFiles(fileNamePrefixdGdz);
1701 }
1702 }
1703
1704 log << "*****************************************"
1705 "**************************************\n";
1706
1707 return;
1708}
Struct to store information on neighbor atoms.
Definition: Atom.h:35
Storage for a single atom.
Definition: Atom.h:32
std::vector< Neighbor > neighbors
Neighbor array (maximum number defined in macros.h.
Definition: Atom.h:148
std::vector< double > G
Symmetry function values.
Definition: Atom.h:134

References combineFiles(), comm, nnp::Mode::elementMap, nnp::Atom::G, nnp::Mode::log, myRank, nnp::Atom::neighbors, nnp::Mode::numElements, nnp::strpr(), and structures.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ collectError()

void Dataset::collectError ( std::string const &  property,
std::map< std::string, double > &  error,
std::size_t &  count 
) const

Collect error metrics of a property over all MPI procs.

Parameters
[in]propertyOne of "energy", "force" or "charge".
[in,out]errorMetric sums of this proc (in), global metric (out).
[in,out]countCount for this proc (in), global count (out).

Definition at line 1710 of file Dataset.cpp.

1713{
1714 if (property == "energy")
1715 {
1716 MPI_Allreduce(MPI_IN_PLACE, &count , 1, MPI_SIZE_T, MPI_SUM, comm);
1717 MPI_Allreduce(MPI_IN_PLACE, &(error.at("RMSEpa")), 1, MPI_DOUBLE, MPI_SUM, comm);
1718 MPI_Allreduce(MPI_IN_PLACE, &(error.at("RMSE")) , 1, MPI_DOUBLE, MPI_SUM, comm);
1719 MPI_Allreduce(MPI_IN_PLACE, &(error.at("MAEpa")) , 1, MPI_DOUBLE, MPI_SUM, comm);
1720 MPI_Allreduce(MPI_IN_PLACE, &(error.at("MAE")) , 1, MPI_DOUBLE, MPI_SUM, comm);
1721 error.at("RMSEpa") = sqrt(error.at("RMSEpa") / count);
1722 error.at("RMSE") = sqrt(error.at("RMSE") / count);
1723 error.at("MAEpa") = error.at("MAEpa") / count;
1724 error.at("MAE") = error.at("MAE") / count;
1725 }
1726 else if (property == "force" || property == "charge")
1727 {
1728 MPI_Allreduce(MPI_IN_PLACE, &count , 1, MPI_SIZE_T, MPI_SUM, comm);
1729 MPI_Allreduce(MPI_IN_PLACE, &(error.at("RMSE")), 1, MPI_DOUBLE, MPI_SUM, comm);
1730 MPI_Allreduce(MPI_IN_PLACE, &(error.at("MAE")) , 1, MPI_DOUBLE, MPI_SUM, comm);
1731 error.at("RMSE") = sqrt(error.at("RMSE") / count);
1732 error.at("MAE") = error.at("MAE") / count;
1733 }
1734 else
1735 {
1736 throw runtime_error("ERROR: Unknown property for error collection.\n");
1737 }
1738
1739 return;
1740}

References comm, and MPI_SIZE_T.

Referenced by nnp::Training::calculateError(), and main().

Here is the caller graph for this function:

◆ combineFiles()

void Dataset::combineFiles ( std::string  filePrefix) const

Combine individual MPI proc files to one.

Parameters
[in]filePrefixFile prefix without the ".0001" suffix.

CAUTION: Make sure files are completely written and closed.

Definition at line 1742 of file Dataset.cpp.

1743{
1744 ofstream combinedFile(filePrefix.c_str(), ios::binary);
1745 if (!combinedFile.is_open())
1746 {
1747 throw runtime_error(strpr("ERROR: Could not open file: %s.\n",
1748 filePrefix.c_str()));
1749 }
1750 for (int i = 0; i < numProcs; ++i)
1751 {
1752 string procFileName = strpr("%s.%04d", filePrefix.c_str(), i);
1753 ifstream procFile(procFileName.c_str(), ios::binary);
1754 if (!procFile.is_open())
1755 {
1756 throw runtime_error(strpr("ERROR: Could not open file: %s.\n",
1757 procFileName.c_str()));
1758 }
1759 // If file is empty, do not get rdbuf, otherwise combined file will be
1760 // closed!
1761 if (procFile.peek() != ifstream::traits_type::eof())
1762 {
1763 combinedFile << procFile.rdbuf();
1764 }
1765 procFile.close();
1766 remove(procFileName.c_str());
1767 }
1768 combinedFile.close();
1769
1770 return;
1771}

References numProcs, and nnp::strpr().

Referenced by nnp::Training::calculateError(), nnp::Training::dataSetNormalization(), main(), writeAtomicEnvironmentFile(), writeNeighborHistogram(), writeNeighborLists(), and nnp::Training::writeSetsToFiles().

Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ structures

◆ myRank

◆ numProcs

int nnp::Dataset::numProcs
protected

◆ numStructures

std::size_t nnp::Dataset::numStructures
protected

Total number of structures in dataset.

Definition at line 203 of file Dataset.h.

Referenced by nnp::Training::dataSetNormalization(), distributeStructures(), prepareNumericForces(), and writeSymmetryFunctionFile().

◆ myName

std::string nnp::Dataset::myName
protected

My processor name.

Definition at line 205 of file Dataset.h.

Referenced by setupMPI().

◆ comm

◆ rng

gsl_rng* nnp::Dataset::rng
protected

GSL random number generator (different seed for each MPI process).

Definition at line 209 of file Dataset.h.

Referenced by distributeStructures(), nnp::Training::selectSets(), setupRandomNumberGenerator(), nnp::Training::setupTraining(), and ~Dataset().

◆ rngGlobal

gsl_rng* nnp::Dataset::rngGlobal
protected

Global GSL random number generator (equal seed for each MPI process).

Definition at line 211 of file Dataset.h.

Referenced by nnp::Training::randomizeNeuralNetworkWeights(), setupRandomNumberGenerator(), nnp::Training::setupTraining(), and ~Dataset().


The documentation for this class was generated from the following files: