modernize electron stopping file reader

This commit is contained in:
Axel Kohlmeyer
2022-01-18 16:48:26 -05:00
parent ef13455d6b
commit b890b564ca

View File

@ -29,6 +29,7 @@
#include "neigh_list.h" #include "neigh_list.h"
#include "neigh_request.h" #include "neigh_request.h"
#include "neighbor.h" #include "neighbor.h"
#include "potential_file_reader.h"
#include "region.h" #include "region.h"
#include "update.h" #include "update.h"
@ -236,49 +237,43 @@ double FixElectronStopping::compute_scalar()
return SeLoss_all; return SeLoss_all;
} }
/* ---------------------------------------------------------------------- */ /* ----------------------------------------------------------------------
read electron stopping parameters. only called from MPI rank 0.
format: energy then one column per atom type
read as many lines as available.
energies must be sorted in ascending order.
---------------------------------------------------------------------- */
void FixElectronStopping::read_table(const char *file) void FixElectronStopping::read_table(const char *file)
{ {
char line[MAXLINE];
FILE *fp = utils::open_potential(file,lmp,nullptr);
if (fp == nullptr)
error->one(FLERR,"Cannot open stopping range table {}: {}", file, utils::getsyserror());
const int ncol = atom->ntypes + 1; const int ncol = atom->ntypes + 1;
int nlines = 0;
PotentialFileReader reader(lmp, file, "electron stopping data table");
int l = 0; try {
while (true) { char *line;
if (fgets(line, MAXLINE, fp) == nullptr) break; // end of file double oldvalue = 0.0;
if (line[0] == '#') continue; // comment
char *pch = strtok(line, " \t\n\r"); while ((line = reader.next_line())) {
if (pch == nullptr) continue; // blank line if (nlines >= maxlines) grow_table();
ValueTokenizer values(line);
elstop_ranges[0][nlines] = values.next_double();
if (elstop_ranges[0][nlines] <= oldvalue)
throw TokenizerException("energy values must be positive and in ascending order",line);
if (l >= maxlines) grow_table(); oldvalue = elstop_ranges[0][nlines];
for (int i = 1; i < ncol; ++i)
elstop_ranges[i][nlines] = values.next_double();
int i = 0; ++nlines;
for ( ; i < ncol && pch != nullptr; i++) {
elstop_ranges[i][l] = utils::numeric(FLERR, pch,false,lmp);
pch = strtok(nullptr, " \t\n\r");
} }
} catch (std::exception &e) {
if (i != ncol || pch != nullptr) // too short or too long error->one(FLERR, "Problem parsing electron stopping data: {}", e.what());
error->one(FLERR, "fix electron/stopping: Invalid table line");
if (l >= 1 && elstop_ranges[0][l] <= elstop_ranges[0][l-1])
error->one(FLERR,
"fix electron/stopping: Energies must be in ascending order");
l++;
} }
table_entries = l; if (nlines == 0)
if (table_entries == 0)
error->one(FLERR, "Did not find any data in electron/stopping table file"); error->one(FLERR, "Did not find any data in electron/stopping table file");
fclose(fp); table_entries = nlines;
} }
/* ---------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- */