git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@8164 f3b2605a-c512-4ea7-a41b-209d697bcdaa

This commit is contained in:
sjplimp
2012-05-23 16:24:30 +00:00
parent 44904c5bf3
commit aa8aeb0294
15 changed files with 109 additions and 23 deletions

View File

@ -71,8 +71,9 @@ support = ["Makefile","Make.sh","Makefile.package.empty",
# packages that have external libs with their external lib dir
extlibs = {"USER-ATC": "atc", "USER-AWPMD": "awpmd", "USER-CUDA": "cuda",
"GPU": "gpu","MEAM": "meam", "POEMS": "poems", "REAX": "reax"}
extlibs = {"USER-ATC": "atc", "USER-AWPMD": "awpmd", "USER-COLVARS": "colvars",
"USER-CUDA": "cuda","GPU": "gpu","MEAM": "meam", "POEMS": "poems",
"REAX": "reax"}
# help messages

View File

@ -17,10 +17,10 @@ PACKAGE = asphere class2 colloid dipole fld gpu granular kim \
kspace manybody mc meam molecule opt peri poems reax replica \
shock srd xtc
PACKUSER = user-misc user-atc user-awpmd user-cg-cmm \
PACKUSER = user-misc user-atc user-awpmd user-cg-cmm user-colvars \
user-cuda user-eff user-ewaldn user-omp user-reaxc user-sph
PACKLIB = gpu kim meam poems reax user-atc user-awpmd user-cuda
PACKLIB = gpu kim meam poems reax user-atc user-awpmd user-colvars user-cuda
PACKALL = $(PACKAGE) $(PACKUSER)

View File

@ -151,6 +151,7 @@ class Fix : protected Pointers {
virtual bigint read_data_skip_lines(char *) {return 0;}
virtual int modify_param(int, char **) {return 0;}
virtual void *extract(const char *, int &) {return NULL;}
virtual double memory_usage() {return 0.0;}

View File

@ -66,6 +66,7 @@ FixLangevin::FixLangevin(LAMMPS *lmp, int narg, char **arg) :
strcpy(tstr,&arg[3][2]);
} else {
t_start = atof(arg[3]);
t_target = t_start;
tstyle = CONSTANT;
}
@ -274,7 +275,7 @@ void FixLangevin::post_force_respa(int vflag, int ilevel, int iloop)
void FixLangevin::post_force_no_tally()
{
double gamma1,gamma2,t_target;
double gamma1,gamma2;
double **v = atom->v;
double **f = atom->f;
@ -461,7 +462,7 @@ void FixLangevin::post_force_no_tally()
void FixLangevin::post_force_tally()
{
double gamma1,gamma2,t_target;
double gamma1,gamma2;
// reallocate flangevin if necessary
@ -727,7 +728,7 @@ void FixLangevin::end_of_step()
void FixLangevin::reset_target(double t_new)
{
t_start = t_stop = t_new;
t_target = t_start = t_stop = t_new;
}
/* ---------------------------------------------------------------------- */
@ -798,6 +799,19 @@ double FixLangevin::compute_scalar()
return -energy_all;
}
/* ----------------------------------------------------------------------
extract thermostat properties
------------------------------------------------------------------------- */
void *FixLangevin::extract(const char *str, int &dim)
{
dim=0;
if (strcmp(str,"t_target") == 0) {
return &t_target;
}
return NULL;
}
/* ----------------------------------------------------------------------
memory usage of tally array
------------------------------------------------------------------------- */

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -39,10 +39,11 @@ class FixLangevin : public Fix {
int modify_param(int, char **);
virtual double compute_scalar();
double memory_usage();
virtual void *extract(const char *, int &);
protected:
int which,tally,zeroflag,oflag,aflag;
double t_start,t_stop,t_period;
double t_start,t_stop,t_period,t_target;
double *gfactor1,*gfactor2,*ratio;
double energy,energy_onestep;
double tsqrt;

View File

@ -115,6 +115,7 @@ FixNH::FixNH(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
if (iarg+4 > narg) error->all(FLERR,"Illegal fix nvt/npt/nph command");
tstat_flag = 1;
t_start = atof(arg[iarg+1]);
t_target = t_start;
t_stop = atof(arg[iarg+2]);
t_period = atof(arg[iarg+3]);
if (t_start < 0.0 || t_stop <= 0.0)
@ -1597,7 +1598,7 @@ double FixNH::compute_vector(int n)
void FixNH::reset_target(double t_new)
{
t_start = t_stop = t_new;
t_target = t_start = t_stop = t_new;
}
/* ---------------------------------------------------------------------- */
@ -1623,6 +1624,19 @@ void FixNH::reset_dt()
tdrag_factor = 1.0 - (update->dt * t_freq * drag / nc_tchain);
}
/* ----------------------------------------------------------------------
extract thermostat properties
------------------------------------------------------------------------- */
void *FixNH::extract(const char *str, int &dim)
{
dim=0;
if (strcmp(str,"t_target") == 0) {
return &t_target;
}
return NULL;
}
/* ----------------------------------------------------------------------
perform half-step update of chain thermostat variables
------------------------------------------------------------------------- */

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -38,6 +38,7 @@ class FixNH : public Fix {
int modify_param(int, char **);
void reset_target(double);
void reset_dt();
virtual void *extract(const char*,int &);
protected:
int dimension,which;

View File

@ -352,6 +352,7 @@ FixRigid::FixRigid(LAMMPS *lmp, int narg, char **arg) :
} else error->all(FLERR,"Illegal fix rigid command");
}
t_target = t_start;
// initialize Marsaglia RNG with processor-unique seed
@ -1231,7 +1232,7 @@ void FixRigid::post_force(int vflag)
double delta = update->ntimestep - update->beginstep;
delta /= update->endstep - update->beginstep;
double t_target = t_start + delta * (t_stop-t_start);
t_target = t_start + delta * (t_stop-t_start);
double tsqrt = sqrt(t_target);
double boltz = force->boltz;
@ -2084,6 +2085,19 @@ double FixRigid::compute_scalar()
return t;
}
/* ----------------------------------------------------------------------
extract thermostat properties
------------------------------------------------------------------------- */
void *FixRigid::extract(const char *str, int &dim)
{
dim=0;
if (strcmp(str,"t_target") == 0) {
return &t_target;
}
return NULL;
}
/* ----------------------------------------------------------------------
return attributes of a rigid body
15 values per body

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -49,6 +49,7 @@ class FixRigid : public Fix {
int dof(int);
void deform(int);
void reset_dt();
virtual void *extract(const char*,int &);
double compute_array(int, int);
protected:
@ -94,7 +95,7 @@ class FixRigid : public Fix {
int langflag; // 0/1 = no/yes Langevin thermostat
int tempflag; // NVT settings
double t_start,t_stop;
double t_start,t_stop,t_target;
double t_period,t_freq;
int t_chain,t_iter,t_order;

View File

@ -44,6 +44,7 @@ FixTempBerendsen::FixTempBerendsen(LAMMPS *lmp, int narg, char **arg) :
extscalar = 1;
t_start = atof(arg[3]);
t_target = t_start;
t_stop = atof(arg[4]);
t_period = atof(arg[5]);
@ -113,7 +114,7 @@ void FixTempBerendsen::end_of_step()
double delta = update->ntimestep - update->beginstep;
delta /= update->endstep - update->beginstep;
double t_target = t_start + delta * (t_stop-t_start);
t_target = t_start + delta * (t_stop-t_start);
// rescale velocities by lamda
// for BIAS:
@ -181,7 +182,7 @@ int FixTempBerendsen::modify_param(int narg, char **arg)
void FixTempBerendsen::reset_target(double t_new)
{
t_start = t_stop = t_new;
t_target = t_start = t_stop = t_new;
}
/* ---------------------------------------------------------------------- */
@ -190,3 +191,17 @@ double FixTempBerendsen::compute_scalar()
{
return energy;
}
/* ----------------------------------------------------------------------
extract thermostat properties
------------------------------------------------------------------------- */
void *FixTempBerendsen::extract(const char *str, int &dim)
{
dim=0;
if (strcmp(str,"t_target") == 0) {
return &t_target;
}
return NULL;
}

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -34,10 +34,11 @@ class FixTempBerendsen : public Fix {
int modify_param(int, char **);
void reset_target(double);
double compute_scalar();
virtual void *extract(const char *, int &);
private:
int which;
double t_start,t_stop,t_period;
double t_start,t_stop,t_period,t_target;
double energy;
char *id_temp;

View File

@ -46,6 +46,7 @@ FixTempRescale::FixTempRescale(LAMMPS *lmp, int narg, char **arg) :
extscalar = 1;
t_start = atof(arg[4]);
t_target = t_start;
t_stop = atof(arg[5]);
t_window = atof(arg[6]);
fraction = atof(arg[7]);
@ -112,7 +113,7 @@ void FixTempRescale::end_of_step()
double delta = update->ntimestep - update->beginstep;
delta /= update->endstep - update->beginstep;
double t_target = t_start + delta * (t_stop-t_start);
t_target = t_start + delta * (t_stop-t_start);
// rescale velocity of appropriate atoms if outside window
// for BIAS:
@ -186,7 +187,7 @@ int FixTempRescale::modify_param(int narg, char **arg)
void FixTempRescale::reset_target(double t_new)
{
t_start = t_stop = t_new;
t_target = t_start = t_stop = t_new;
}
/* ---------------------------------------------------------------------- */
@ -195,3 +196,17 @@ double FixTempRescale::compute_scalar()
{
return energy;
}
/* ----------------------------------------------------------------------
extract thermostat properties
------------------------------------------------------------------------- */
void *FixTempRescale::extract(const char *str, int &dim)
{
dim=0;
if (strcmp(str,"t_target") == 0) {
return &t_target;
}
return NULL;
}

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -34,10 +34,11 @@ class FixTempRescale : public Fix {
int modify_param(int, char **);
void reset_target(double);
double compute_scalar();
virtual void *extract(const char *, int &);
protected:
int which;
double t_start,t_stop,t_window;
double t_start,t_stop,t_window,t_target;
double fraction,energy,efactor;
char *id_temp;

View File

@ -1,4 +1,4 @@
/* ----------------------------------------------------------------------
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
@ -37,6 +37,7 @@ class Force : protected Pointers {
double mvh2r; // conversion of mv/hbar to distance
// hbar = h/(2*pi)
double angstrom; // 1 angstrom in native units
double femtosecond; // 1 femtosecond in native units
double qelectron; // 1 electron charge abs() in native units
int newton,newton_pair,newton_bond; // Newton's 3rd law settings

View File

@ -132,6 +132,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0;
force->mvh2r = 0.0;
force->angstrom = 1.0;
force->femtosecond = 1.0;
force->qelectron = 1.0;
dt = 0.005;
@ -152,6 +153,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0957018663603261;
force->mvh2r = 1.5339009481951;
force->angstrom = 1.0;
force->femtosecond = 1.0;
force->qelectron = 1.0;
dt = 1.0;
@ -172,6 +174,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0;
force->mvh2r = 0.0;
force->angstrom = 1.0;
force->femtosecond = 1.0e-3;
force->qelectron = 1.0;
dt = 0.001;
@ -192,6 +195,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0;
force->mvh2r = 0.0;
force->angstrom = 1.0e-10;
force->femtosecond = 1.0e-15;
force->qelectron = 1.6021765e-19;
dt = 1.0e-8;
@ -212,6 +216,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0;
force->mvh2r = 0.0;
force->angstrom = 1.0e-8;
force->femtosecond = 1.0e-15;
force->qelectron = 4.8032044e-10;
dt = 1.0e-8;
@ -232,6 +237,7 @@ void Update::set_units(const char *style)
force->hhmrr2e = 0.0;
force->mvh2r = 0.0;
force->angstrom = 1.88972612;
force->femtosecond = 0.0241888428;
force->qelectron = 1.0;
dt = 0.001;