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

This commit is contained in:
sjplimp
2010-01-13 20:55:06 +00:00
parent 363e6c371f
commit 1dbad849b7
7 changed files with 39 additions and 29 deletions

View File

@ -52,6 +52,10 @@ Atom::Atom(LAMMPS *lmp) : Pointers(lmp)
extra_bond_per_atom = 0; extra_bond_per_atom = 0;
firstgroupname = NULL; firstgroupname = NULL;
sortfreq = 0;
maxbin = maxnext = 0;
binhead = NULL;
next = permute = NULL;
// initialize atom arrays // initialize atom arrays
// customize by adding new array // customize by adding new array
@ -110,13 +114,6 @@ Atom::Atom(LAMMPS *lmp) : Pointers(lmp)
nextra_store = 0; nextra_store = 0;
extra = NULL; extra = NULL;
// sorting
sortflag = 0;
maxbin = maxnext = 0;
binhead = NULL;
next = permute = NULL;
// default mapping values and hash table primes // default mapping values and hash table primes
tag_enable = 1; tag_enable = 1;
@ -146,7 +143,11 @@ Atom::~Atom()
{ {
delete [] atom_style; delete [] atom_style;
delete avec; delete avec;
delete [] firstgroupname; delete [] firstgroupname;
memory->sfree(binhead);
memory->sfree(next);
memory->sfree(permute);
// delete atom arrays // delete atom arrays
// customize by adding new array // customize by adding new array
@ -211,12 +212,6 @@ Atom::~Atom()
delete [] dipole; delete [] dipole;
delete [] dipole_setflag; delete [] dipole_setflag;
// delete sorting arrays
memory->sfree(binhead);
memory->sfree(next);
memory->sfree(permute);
// delete extra arrays // delete extra arrays
memory->sfree(extra_grow); memory->sfree(extra_grow);
@ -315,7 +310,7 @@ void Atom::setup()
// setup for sorting // setup for sorting
// binsize = user setting or 1/2 of neighbor cutoff // binsize = user setting or 1/2 of neighbor cutoff
if (sortflag) { if (sortfreq > 0) {
double binsize; double binsize;
if (userbinsize > 0.0) binsize = userbinsize; if (userbinsize > 0.0) binsize = userbinsize;
else binsize = 0.5 * neighbor->cutneighmax; else binsize = 0.5 * neighbor->cutneighmax;
@ -341,6 +336,8 @@ int Atom::style_match(const char *style)
/* ---------------------------------------------------------------------- /* ----------------------------------------------------------------------
modify parameters of the atom style modify parameters of the atom style
some options can only be invoked before simulation box is defined
first and sort options cannot be used together
------------------------------------------------------------------------- */ ------------------------------------------------------------------------- */
void Atom::modify_params(int narg, char **arg) void Atom::modify_params(int narg, char **arg)
@ -354,6 +351,8 @@ void Atom::modify_params(int narg, char **arg)
if (strcmp(arg[iarg+1],"array") == 0) map_style = 1; if (strcmp(arg[iarg+1],"array") == 0) map_style = 1;
else if (strcmp(arg[iarg+1],"hash") == 0) map_style = 2; else if (strcmp(arg[iarg+1],"hash") == 0) map_style = 2;
else error->all("Illegal atom_modify command"); else error->all("Illegal atom_modify command");
if (domain->box_exist)
error->all("Atom_modify map command after simulation box is defined");
iarg += 2; iarg += 2;
} else if (strcmp(arg[iarg],"first") == 0) { } else if (strcmp(arg[iarg],"first") == 0) {
if (iarg+2 > narg) error->all("Illegal atom_modify command"); if (iarg+2 > narg) error->all("Illegal atom_modify command");
@ -363,6 +362,7 @@ void Atom::modify_params(int narg, char **arg)
firstgroupname = new char[n]; firstgroupname = new char[n];
strcpy(firstgroupname,arg[iarg+1]); strcpy(firstgroupname,arg[iarg+1]);
} }
sortfreq = 0;
iarg += 2; iarg += 2;
} else if (strcmp(arg[iarg],"sort") == 0) { } else if (strcmp(arg[iarg],"sort") == 0) {
if (iarg+3 > narg) error->all("Illegal atom_modify command"); if (iarg+3 > narg) error->all("Illegal atom_modify command");
@ -370,8 +370,9 @@ void Atom::modify_params(int narg, char **arg)
userbinsize = atof(arg[iarg+2]); userbinsize = atof(arg[iarg+2]);
if (sortfreq < 0 || userbinsize < 0.0) if (sortfreq < 0 || userbinsize < 0.0)
error->all("Illegal atom_modify command"); error->all("Illegal atom_modify command");
if (sortfreq == 0) sortflag = 0; if (sortfreq >= 0 && firstgroupname)
else sortflag = 1; error->all("Atom_modify sort and first options "
"cannot be used together");
iarg += 3; iarg += 3;
} else error->all("Illegal atom_modify command"); } else error->all("Illegal atom_modify command");
} }
@ -1370,7 +1371,7 @@ void Atom::sort()
permute = (int *) memory->smalloc(maxnext*sizeof(int),"atom:permute"); permute = (int *) memory->smalloc(maxnext*sizeof(int),"atom:permute");
} }
// insure one extra location at end of atom arrays // insure there is one extra atom location at end of arrays for swaps
if (nlocal == nmax) avec->grow(0); if (nlocal == nmax) avec->grow(0);

View File

@ -102,7 +102,7 @@ class Atom : protected Pointers {
// spatial sorting of atoms // spatial sorting of atoms
int sortflag; // 0 = off, 1 = on int sortfreq; // sort atoms every this many steps, 0 = off
int nextsort; // next timestep to sort on int nextsort; // next timestep to sort on
// functions // functions
@ -203,7 +203,6 @@ class Atom : protected Pointers {
int *permute; // permutation vector int *permute; // permutation vector
double userbinsize; // sorting bin size double userbinsize; // sorting bin size
double bininv; double bininv;
int sortfreq;
int memlength; // allocated size of memstr int memlength; // allocated size of memstr
char *memstr; // string of array names already counted char *memstr; // string of array names already counted

View File

@ -172,7 +172,7 @@ double ComputePressure::compute_scalar()
double t; double t;
if (keflag) { if (keflag) {
if (temperature->invoked_scalar == update->ntimestep) if (temperature->invoked_scalar != update->ntimestep)
t = temperature->compute_scalar(); t = temperature->compute_scalar();
else t = temperature->scalar; else t = temperature->scalar;
} }
@ -209,7 +209,7 @@ void ComputePressure::compute_vector()
if (update->vflag_global != invoked_vector) if (update->vflag_global != invoked_vector)
error->all("Virial was not tallied on needed timestep"); error->all("Virial was not tallied on needed timestep");
// invoke temperature it it hasn't been already // invoke temperature if it hasn't been already
double *ke_tensor; double *ke_tensor;
if (keflag) { if (keflag) {

View File

@ -691,8 +691,6 @@ void Input::angle_style()
void Input::atom_modify() void Input::atom_modify()
{ {
if (domain->box_exist)
error->all("Atom_modify command after simulation box is defined");
atom->modify_params(narg,arg); atom->modify_params(narg,arg);
} }

View File

@ -192,12 +192,20 @@ void Modify::init()
list_init_compute(); list_init_compute();
// init each compute // init each compute
// set invoked_scalar,vector,etc to -1 to force new run to re-compute them
// add initial timestep to all computes that store invocation times // add initial timestep to all computes that store invocation times
// since any of them may be invoked by initial thermo // since any of them may be invoked by initial thermo
// do not clear out invocation times stored within a compute, // do not clear out invocation times stored within a compute,
// b/c some may be holdovers from previous run, like for ave fixes // b/c some may be holdovers from previous run, like for ave fixes
for (i = 0; i < ncompute; i++) compute[i]->init(); for (i = 0; i < ncompute; i++) {
compute[i]->init();
compute[i]->invoked_scalar = -1;
compute[i]->invoked_vector = -1;
compute[i]->invoked_array = -1;
compute[i]->invoked_peratom = -1;
compute[i]->invoked_local = -1;
}
addstep_compute_all(update->ntimestep); addstep_compute_all(update->ntimestep);
// warn if any particle is time integrated more than once // warn if any particle is time integrated more than once

View File

@ -222,7 +222,8 @@ void Update::create_minimize(int narg, char **arg)
do not allow any timestep-dependent fixes to be defined do not allow any timestep-dependent fixes to be defined
do not allow any dynamic regions to be defined do not allow any dynamic regions to be defined
reset eflag/vflag global so nothing will think eng/virial are current reset eflag/vflag global so nothing will think eng/virial are current
reset invoked flags of computes, so nothing will think they are current reset invoked flags of computes,
so nothing will think they are current between runs
clear timestep list of computes that store future invocation times clear timestep list of computes that store future invocation times
------------------------------------------------------------------------- */ ------------------------------------------------------------------------- */
@ -249,6 +250,7 @@ void Update::reset_timestep(int narg, char **arg)
for (int i = 0; i < modify->ncompute; i++) { for (int i = 0; i < modify->ncompute; i++) {
modify->compute[i]->invoked_scalar = -1; modify->compute[i]->invoked_scalar = -1;
modify->compute[i]->invoked_vector = -1; modify->compute[i]->invoked_vector = -1;
modify->compute[i]->invoked_array = -1;
modify->compute[i]->invoked_peratom = -1; modify->compute[i]->invoked_peratom = -1;
modify->compute[i]->invoked_local = -1; modify->compute[i]->invoked_local = -1;
} }

View File

@ -92,7 +92,7 @@ void Verlet::setup()
comm->setup(); comm->setup();
if (neighbor->style) neighbor->setup_bins(); if (neighbor->style) neighbor->setup_bins();
comm->exchange(); comm->exchange();
if (atom->sortflag) atom->sort(); if (atom->sortfreq > 0) atom->sort();
comm->borders(); comm->borders();
if (triclinic) domain->lamda2x(atom->nlocal+atom->nghost); if (triclinic) domain->lamda2x(atom->nlocal+atom->nghost);
neighbor->build(); neighbor->build();
@ -178,7 +178,7 @@ void Verlet::setup_minimal(int flag)
void Verlet::run(int n) void Verlet::run(int n)
{ {
int nflag,ntimestep; int nflag,ntimestep,sortflag;
int n_post_integrate = modify->n_post_integrate; int n_post_integrate = modify->n_post_integrate;
int n_pre_exchange = modify->n_pre_exchange; int n_pre_exchange = modify->n_pre_exchange;
@ -186,7 +186,9 @@ void Verlet::run(int n)
int n_pre_force = modify->n_pre_force; int n_pre_force = modify->n_pre_force;
int n_post_force = modify->n_post_force; int n_post_force = modify->n_post_force;
int n_end_of_step = modify->n_end_of_step; int n_end_of_step = modify->n_end_of_step;
int sortflag = atom->sortflag;
if (atom->sortfreq > 0) sortflag = 1;
else sortflag = 0;
for (int i = 0; i < n; i++) { for (int i = 0; i < n; i++) {
@ -217,7 +219,7 @@ void Verlet::run(int n)
} }
timer->stamp(); timer->stamp();
comm->exchange(); comm->exchange();
if (sortflag && ntimestep > atom->nextsort) atom->sort(); if (sortflag && ntimestep >= atom->nextsort) atom->sort();
comm->borders(); comm->borders();
if (triclinic) domain->lamda2x(atom->nlocal+atom->nghost); if (triclinic) domain->lamda2x(atom->nlocal+atom->nghost);
timer->stamp(TIME_COMM); timer->stamp(TIME_COMM);