support offsets for molecule IDs (if available) in read_data similar to atomIDs

suggested by felipe perez in https://sourceforge.net/p/lammps/mailman/message/36236631/
This commit is contained in:
Axel Kohlmeyer
2018-02-23 18:02:05 -05:00
parent 0003bb6766
commit bba4bd1489
5 changed files with 47 additions and 25 deletions

View File

@ -15,10 +15,11 @@ read_data file keyword args ... :pre
file = name of data file to read in :ulb,l
zero or more keyword/arg pairs may be appended :l
keyword = {add} or {offset} or {shift} or {extra/atom/types} or {extra/bond/types} or {extra/angle/types} or {extra/dihedral/types} or {extra/improper/types} or {extra/bond/per/atom} or {extra/angle/per/atom} or {extra/dihedral/per/atom} or {extra/improper/per/atom} or {group} or {nocoeff} or {fix} :l
{add} arg = {append} or {Nstart} or {merge}
append = add new atoms with IDs appended to current IDs
Nstart = add new atoms with IDs starting with Nstart
merge = add new atoms with their IDs unchanged
{add} arg = {append} or {IDoffset} or {IDoffset MOLoffset} or {merge}
append = add new atoms with atom IDs appended to current IDs
IDoffset = add new atoms with atom IDs having IDoffset added
MOLoffset = add new atoms with molecule IDs having MOLoffset added (only when molecule IDs are enabled)
merge = add new atoms with their atom IDs (and molecule IDs) unchanged
{offset} args = toff boff aoff doff ioff
toff = offset to add to atom types
boff = offset to add to bond types
@ -120,20 +121,26 @@ boundary, then the atoms may become far apart if the box size grows.
This will separate the atoms in the bond, which can lead to "lost"
bond atoms or bad dynamics.
The three choices for the {add} argument affect how the IDs of atoms
in the data file are treated. If {append} is specified, atoms in the
data file are added to the current system, with their atom IDs reset
so that an atomID = M in the data file becomes atomID = N+M, where N
is the largest atom ID in the current system. This rule is applied to
all occurrences of atom IDs in the data file, e.g. in the Velocity or
Bonds section. If {Nstart} is specified, then {Nstart} is a numeric
value is given, e.g. 1000, so that an atomID = M in the data file
becomes atomID = 1000+M. If {merge} is specified, the data file atoms
The three choices for the {add} argument affect how the atom IDs and
molecule IDs of atoms in the data file are treated. If {append} is
specified, atoms in the data file are added to the current system,
with their atom IDs reset so that an atomID = M in the data file
becomes atomID = N+M, where N is the largest atom ID in the current
system. This rule is applied to all occurrences of atom IDs in the
data file, e.g. in the Velocity or Bonds section. This is also done
for molecule IDs, if the atom style does support molecule IDs or
they are enabled via fix property/atom. If {IDoffset} is specified,
then {IDoffset} is a numeric value is given, e.g. 1000, so that an
atomID = M in the data file becomes atomID = 1000+M. For systems
with enabled molecule IDs, another numerical argument {MOLoffset}
is required representing the equivalent offset for molecule IDs.
If {merge} is specified, the data file atoms
are added to the current system without changing their IDs. They are
assumed to merge (without duplication) with the currently defined
atoms. It is up to you to insure there are no multiply defined atom
IDs, as LAMMPS only performs an incomplete check that this is the case
by insuring the resulting max atomID >= the number of atoms.
by insuring the resulting max atomID >= the number of atoms. For
molecule IDs, there is no check done at all.
The {offset} and {shift} keywords can only be used if the {add}
keyword is also specified.

View File

@ -822,8 +822,8 @@ void Atom::deallocate_topology()
call style-specific routine to parse line
------------------------------------------------------------------------- */
void Atom::data_atoms(int n, char *buf, tagint id_offset, int type_offset,
int shiftflag, double *shift)
void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
int type_offset, int shiftflag, double *shift)
{
int m,xptr,iptr;
imageint imagedata;
@ -948,6 +948,7 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, int type_offset,
coord[2] >= sublo[2] && coord[2] < subhi[2]) {
avec->data_atom(xdata,imagedata,values);
if (id_offset) tag[nlocal-1] += id_offset;
if (mol_offset) molecule[nlocal-1] += mol_offset;
if (type_offset) {
type[nlocal-1] += type_offset;
if (type[nlocal-1] > ntypes)

View File

@ -229,7 +229,7 @@ class Atom : protected Pointers {
void deallocate_topology();
void data_atoms(int, char *, tagint, int, int, double *);
void data_atoms(int, char *, tagint, tagint, int, int, double *);
void data_vels(int, char *, tagint);
void data_bonds(int, char *, int *, tagint, int);
void data_angles(int, char *, int *, tagint, int);

View File

@ -124,7 +124,7 @@ void ReadData::command(int narg, char **arg)
addflag = NONE;
coeffflag = 1;
id_offset = 0;
id_offset = mol_offset = 0;
offsetflag = shiftflag = 0;
toffset = boffset = aoffset = doffset = ioffset = 0;
shift[0] = shift[1] = shift[2] = 0.0;
@ -145,11 +145,21 @@ void ReadData::command(int narg, char **arg)
if (strcmp(arg[iarg+1],"append") == 0) addflag = APPEND;
else if (strcmp(arg[iarg+1],"merge") == 0) addflag = MERGE;
else {
if (atom->molecule_flag && (iarg+3 > narg))
error->all(FLERR,"Illegal read_data command");
addflag = VALUE;
bigint offset = force->bnumeric(FLERR,arg[iarg+1]);
if (offset > MAXTAGINT)
error->all(FLERR,"Read data add offset is too big");
error->all(FLERR,"Read data add atomID offset is too big");
id_offset = offset;
if (atom->molecule_flag) {
offset = force->bnumeric(FLERR,arg[iarg+2]);
if (offset > MAXTAGINT)
error->all(FLERR,"Read data add molID offset is too big");
mol_offset = offset;
iarg++;
}
}
iarg += 2;
} else if (strcmp(arg[iarg],"offset") == 0) {
@ -310,14 +320,18 @@ void ReadData::command(int narg, char **arg)
update->ntimestep = 0;
}
// compute atomID offset for addflag = MERGE
// compute atomID and optionally moleculeID offset for addflag = APPEND
if (addflag == APPEND) {
tagint *tag = atom->tag;
tagint *molecule = atom->molecule;
int nlocal = atom->nlocal;
tagint max = 0;
for (int i = 0; i < nlocal; i++) max = MAX(max,tag[i]);
MPI_Allreduce(&max,&id_offset,1,MPI_LMP_TAGINT,MPI_MAX,world);
tagint maxid = 0, maxmol = 0;
for (int i = 0; i < nlocal; i++) maxid = MAX(maxid,tag[i]);
if (atom->molecule_flag)
for (int i = 0; i < nlocal; i++) maxmol = MAX(maxmol,molecule[i]);
MPI_Allreduce(&maxid,&id_offset,1,MPI_LMP_TAGINT,MPI_MAX,world);
MPI_Allreduce(&maxmol,&mol_offset,1,MPI_LMP_TAGINT,MPI_MAX,world);
}
// set up pointer to hold original styles while we replace them with "zero"
@ -1137,7 +1151,7 @@ void ReadData::atoms()
nchunk = MIN(natoms-nread,CHUNK);
eof = comm->read_lines_from_file(fp,nchunk,MAXLINE,buffer);
if (eof) error->all(FLERR,"Unexpected end of data file");
atom->data_atoms(nchunk,buffer,id_offset,toffset,shiftflag,shift);
atom->data_atoms(nchunk,buffer,id_offset,mol_offset,toffset,shiftflag,shift);
nread += nchunk;
}

View File

@ -39,7 +39,7 @@ class ReadData : protected Pointers {
int narg,maxarg;
char argoffset1[8],argoffset2[8];
bigint id_offset;
bigint id_offset, mol_offset;
int nlocal_previous;
bigint natoms;