coding style updates

This commit is contained in:
Axel Kohlmeyer
2023-12-14 11:15:28 -05:00
parent 727c3c9572
commit 298ce1863c

View File

@ -1,4 +1,3 @@
// clang-format off
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
https://www.lammps.org/, Sandia National Laboratories
@ -31,49 +30,49 @@
using namespace LAMMPS_NS;
using namespace FixConst;
enum{BONDMAX,TLIMIT,DISKFREE,VARIABLE};
enum{LT,LE,GT,GE,EQ,NEQ,XOR};
enum{HARD,SOFT,CONTINUE};
enum{NOMSG=0,YESMSG=1};
enum { BONDMAX, TLIMIT, DISKFREE, VARIABLE };
enum { LT, LE, GT, GE, EQ, NEQ, XOR };
enum { HARD, SOFT, CONTINUE };
enum { NOMSG = 0, YESMSG = 1 };
/* ---------------------------------------------------------------------- */
FixHalt::FixHalt(LAMMPS *lmp, int narg, char **arg) :
Fix(lmp, narg, arg), idvar(nullptr), dlimit_path(nullptr)
Fix(lmp, narg, arg), idvar(nullptr), dlimit_path(nullptr)
{
if (narg < 7) error->all(FLERR,"Illegal fix halt command");
nevery = utils::inumeric(FLERR,arg[3],false,lmp);
if (nevery <= 0) error->all(FLERR,"Illegal fix halt command");
if (narg < 7) utils::missing_cmd_args(FLERR, "fix halt", error);
nevery = utils::inumeric(FLERR, arg[3], false, lmp);
if (nevery <= 0) error->all(FLERR, "Illegal fix halt command: nevery must be > 0");
// comparison args
idvar = nullptr;
int iarg = 4;
if (strcmp(arg[iarg],"tlimit") == 0) {
if (strcmp(arg[iarg], "tlimit") == 0) {
attribute = TLIMIT;
} else if (strcmp(arg[iarg],"diskfree") == 0) {
} else if (strcmp(arg[iarg], "diskfree") == 0) {
attribute = DISKFREE;
dlimit_path = utils::strdup(".");
} else if (strcmp(arg[iarg],"bondmax") == 0) {
} else if (strcmp(arg[iarg], "bondmax") == 0) {
attribute = BONDMAX;
} else {
ArgInfo argi(arg[iarg],ArgInfo::VARIABLE);
ArgInfo argi(arg[iarg], ArgInfo::VARIABLE);
if ((argi.get_type() == ArgInfo::UNKNOWN)
|| (argi.get_type() == ArgInfo::NONE)
|| (argi.get_dim() != 0))
error->all(FLERR,"Invalid fix halt attribute");
if ((argi.get_type() == ArgInfo::UNKNOWN) || (argi.get_type() == ArgInfo::NONE) ||
(argi.get_dim() != 0))
error->all(FLERR, "Invalid fix halt attribute {}", arg[iarg]);
attribute = VARIABLE;
idvar = argi.copy_name();
ivar = input->variable->find(idvar);
if (ivar < 0) error->all(FLERR,"Could not find fix halt variable name");
if (ivar < 0) error->all(FLERR, "Could not find fix halt variable name");
if (input->variable->equalstyle(ivar) == 0)
error->all(FLERR,"Fix halt variable is not equal-style variable");
error->all(FLERR, "Fix halt variable is not equal-style variable");
}
// clang-format off
++iarg;
if (strcmp(arg[iarg],"<") == 0) operation = LT;
else if (strcmp(arg[iarg],"<=") == 0) operation = LE;
@ -85,7 +84,7 @@ FixHalt::FixHalt(LAMMPS *lmp, int narg, char **arg) :
else error->all(FLERR,"Invalid fix halt operator");
++iarg;
value = utils::numeric(FLERR,arg[iarg],false,lmp);
value = utils::numeric(FLERR, arg[iarg], false, lmp);
// parse optional args
@ -93,38 +92,40 @@ FixHalt::FixHalt(LAMMPS *lmp, int narg, char **arg) :
msgflag = YESMSG;
++iarg;
while (iarg < narg) {
if (strcmp(arg[iarg],"error") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix halt command");
if (strcmp(arg[iarg+1],"hard") == 0) eflag = HARD;
else if (strcmp(arg[iarg+1],"soft") == 0) eflag = SOFT;
else if (strcmp(arg[iarg+1],"continue") == 0) eflag = CONTINUE;
else error->all(FLERR,"Illegal fix halt command");
if (strcmp(arg[iarg], "error") == 0) {
if (iarg + 2 > narg) utils::missing_cmd_args(FLERR, "fix halt error", error);
if (strcmp(arg[iarg + 1], "hard") == 0) eflag = HARD;
else if (strcmp(arg[iarg + 1], "soft") == 0) eflag = SOFT;
else if (strcmp(arg[iarg + 1], "continue") == 0) eflag = CONTINUE;
else error->all(FLERR, "Unknown fix halt error condition {}", arg[iarg]);
iarg += 2;
} else if (strcmp(arg[iarg],"message") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix halt command");
msgflag = utils::logical(FLERR,arg[iarg+1],false,lmp);
} else if (strcmp(arg[iarg], "message") == 0) {
if (iarg + 2 > narg) utils::missing_cmd_args(FLERR, "fix halt message", error);
msgflag = utils::logical(FLERR, arg[iarg + 1], false, lmp);
iarg += 2;
} else if (strcmp(arg[iarg],"path") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix halt command");
} else if (strcmp(arg[iarg], "path") == 0) {
if (iarg + 2 > narg) utils::missing_cmd_args(FLERR, "fix halt error", error);
++iarg;
delete[] dlimit_path;
// strip off outer quotes, if present
int len = strlen(arg[iarg])+1;
if ( ((arg[iarg][0] == '"') || (arg[iarg][0] == '\''))
&& (arg[iarg][0] == arg[iarg][len-2])) {
arg[iarg][len-2] = '\0';
dlimit_path = utils::strdup(arg[iarg]+1);
} else dlimit_path = utils::strdup(arg[iarg]);
int len = strlen(arg[iarg]) + 1;
if (((arg[iarg][0] == '"') || (arg[iarg][0] == '\'')) &&
(arg[iarg][0] == arg[iarg][len - 2])) {
arg[iarg][len - 2] = '\0';
dlimit_path = utils::strdup(arg[iarg] + 1);
} else
dlimit_path = utils::strdup(arg[iarg]);
++iarg;
} else error->all(FLERR,"Illegal fix halt command");
} else error->all(FLERR, "Unknown fix halt keyword {}", arg[iarg]);
}
// clang-format on
// add nfirst to all computes that store invocation times
// since don't know a priori which are invoked via variables by this fix
// once in end_of_step() can set timestep for ones actually invoked
if (attribute == VARIABLE) {
const bigint nfirst = (update->ntimestep/nevery)*nevery + nevery;
const bigint nfirst = (update->ntimestep / nevery) * nevery + nevery;
modify->addstep_compute_all(nfirst);
}
}
@ -133,8 +134,8 @@ FixHalt::FixHalt(LAMMPS *lmp, int narg, char **arg) :
FixHalt::~FixHalt()
{
delete [] idvar;
delete [] dlimit_path;
delete[] idvar;
delete[] dlimit_path;
}
/* ---------------------------------------------------------------------- */
@ -156,14 +157,14 @@ void FixHalt::init()
if (attribute == VARIABLE) {
ivar = input->variable->find(idvar);
if (ivar < 0) error->all(FLERR,"Could not find fix halt variable name");
if (ivar < 0) error->all(FLERR, "Could not find fix halt variable {}", idvar);
if (input->variable->equalstyle(ivar) == 0)
error->all(FLERR,"Fix halt variable is not equal-style variable");
error->all(FLERR, "Fix halt variable {} is not equal-style variable", idvar);
}
// settings used by TLIMIT
nextstep = (update->ntimestep/nevery)*nevery + nevery;
nextstep = (update->ntimestep / nevery) * nevery + nevery;
thisstep = -1;
tratio = 0.5;
@ -205,6 +206,10 @@ void FixHalt::end_of_step()
modify->addstep_compute(update->ntimestep + nevery);
}
// ensure that the attribute is *exactly* the same on all ranks
MPI_Bcast(&attvalue, 1, MPI_DOUBLE, 0, world);
// check if halt is triggered, else just return
if (operation == LT) {
@ -220,21 +225,19 @@ void FixHalt::end_of_step()
} else if (operation == NEQ) {
if (attvalue == value) return;
} else if (operation == XOR) {
if ((attvalue == 0.0 && value == 0.0) ||
(attvalue != 0.0 && value != 0.0)) return;
if ((attvalue == 0.0 && value == 0.0) || (attvalue != 0.0 && value != 0.0)) return;
}
// hard halt -> exit LAMMPS
// soft/continue halt -> trigger timer to break from run loop
// print message with ID of fix halt in case multiple instances
std::string message = fmt::format("Fix halt condition for fix-id {} met on "
"step {} with value {}",
std::string message = fmt::format("Fix halt condition for fix-id {} met on step {} with value {}",
id, update->ntimestep, attvalue);
if (eflag == HARD) {
error->all(FLERR,message);
} else if (eflag == SOFT || eflag == CONTINUE) {
if (comm->me == 0 && msgflag == YESMSG) error->message(FLERR,message);
error->all(FLERR, message);
} else if ((eflag == SOFT) || (eflag == CONTINUE)) {
if ((comm->me == 0) && (msgflag == YESMSG)) error->message(FLERR, message);
timer->force_timeout();
}
}
@ -260,8 +263,8 @@ double FixHalt::bondmax()
int **bondlist = neighbor->bondlist;
int nbondlist = neighbor->nbondlist;
int i1,i2;
double delx,dely,delz,rsq;
int i1, i2;
double delx, dely, delz, rsq;
double maxone = 0.0;
for (int n = 0; n < nbondlist; n++) {
@ -272,12 +275,12 @@ double FixHalt::bondmax()
dely = x[i1][1] - x[i2][1];
delz = x[i1][2] - x[i2][2];
rsq = delx*delx + dely*dely + delz*delz;
maxone = MAX(rsq,maxone);
rsq = delx * delx + dely * dely + delz * delz;
maxone = MAX(rsq, maxone);
}
double maxall;
MPI_Allreduce(&maxone,&maxall,1,MPI_DOUBLE,MPI_MAX,world);
MPI_Allreduce(&maxone, &maxall, 1, MPI_DOUBLE, MPI_MAX, world);
return sqrt(maxall);
}
@ -291,48 +294,15 @@ double FixHalt::bondmax()
double FixHalt::tlimit()
{
double cpu = timer->elapsed(Timer::TOTAL);
MPI_Bcast(&cpu,1,MPI_DOUBLE,0,world);
MPI_Bcast(&cpu, 1, MPI_DOUBLE, 0, world);
if (cpu < value) {
bigint elapsed = update->ntimestep - update->firststep;
bigint final = update->firststep +
static_cast<bigint> (tratio*value/cpu * elapsed);
nextstep = (final/nevery)*nevery + nevery;
bigint final = update->firststep + static_cast<bigint>(tratio * value / cpu * elapsed);
nextstep = (final / nevery) * nevery + nevery;
if (nextstep == update->ntimestep) nextstep += nevery;
tratio = 1.0;
}
return cpu;
}
/* ----------------------------------------------------------------------
determine available disk space, if supported. Return -1 if not.
------------------------------------------------------------------------- */
#if defined(__linux__) || defined(__APPLE__) || defined(__FreeBSD__) || defined(__DragonFly__) || defined(__OpenBSD__) || defined(__NetBSD__)
#include <sys/statvfs.h>
#endif
double FixHalt::diskfree()
{
#if defined(__linux__) || defined(__APPLE__) || defined(__FreeBSD__) || defined(__DragonFly__) || defined(__OpenBSD__) || defined(__NetBSD__)
struct statvfs fs;
double disk_free = -1.0;
if (dlimit_path) {
disk_free = 1.0e100;
int rv = statvfs(dlimit_path,&fs);
if (rv == 0) {
#if defined(__linux__)
disk_free = fs.f_bavail*fs.f_bsize/1048576.0;
#elif defined(__APPLE__) || defined(__FreeBSD__) || defined(__DragonFly__) || defined(__OpenBSD__) || defined(__NetBSD__)
disk_free = fs.f_bavail*fs.f_frsize/1048576.0;
#endif
} else
disk_free = -1.0;
MPI_Bcast(&disk_free,1,MPI_DOUBLE,0,world);
}
return disk_free;
#else
return -1.0;
#endif
}