more changes to support late initialization of masses

This commit is contained in:
Axel Kohlmeyer
2013-06-19 00:22:43 +02:00
parent 5a7dd0e74e
commit 7d444dcb71
7 changed files with 40 additions and 39 deletions

View File

@ -669,9 +669,12 @@ void colvar::disable (colvar::task const &t)
void colvar::setup() {
// loop over all components to reset masses of all groups
for (std::vector<cvc *>::iterator cci = cvcs.begin();
cci != cvcs.end(); cci++) {
(*cci)->setup();
for (size_t i = 0; i < cvcs.size(); i++) {
for (size_t ig = 0; ig < cvcs[i]->atom_groups.size(); ig++) {
cvm::atom_group &atoms = *(cvcs[i]->atom_groups[ig]);
atoms.read_positions();
atoms.reset_mass();
}
}
}

View File

@ -55,16 +55,6 @@ void colvar::cvc::parse_group (std::string const &conf,
}
void colvar::cvc::setup()
{
// loop over all atom groups to reset their masses
for (std::vector<cvm::atom_group *>::iterator agi = atom_groups.begin();
agi != atom_groups.end(); agi++) {
(*agi)->reset_mass();
}
}
colvar::cvc::~cvc()
{}

View File

@ -113,10 +113,6 @@ public:
/// objects are declared within other ones)
cvc();
/// \brief Reinitialize internal data for MD codes that
/// can change parameters between runs
void setup();
/// Destructor
virtual ~cvc();

View File

@ -2,7 +2,7 @@
#define COLVARMODULE_H
#ifndef COLVARS_VERSION
#define COLVARS_VERSION "2013-06-18"
#define COLVARS_VERSION "2013-06-19"
#endif
#ifndef COLVARS_DEBUG

View File

@ -97,12 +97,11 @@ colvarproxy_lammps::colvarproxy_lammps(LAMMPS_NS::LAMMPS *lmp,
}
void colvarproxy_lammps::init(const char *conf_file, const int *typemap)
void colvarproxy_lammps::init(const char *conf_file)
{
_typemap = typemap;
// create the colvarmodule instance
colvars = new colvarmodule(conf_file,this);
cvm::log ("Notice! Postponing initializing of masses.\n");
if (_lmp->update->ntimestep != 0) {
cvm::log ("Initializing step number as firstTimestep.\n");
@ -117,9 +116,6 @@ void colvarproxy_lammps::init(const char *conf_file, const int *typemap)
log(cvm::line_marker);
log("Info: done initializing the colvars proxy object.\n");
}
// this is only valid in this function.
_typemap = NULL;
}
colvarproxy_lammps::~colvarproxy_lammps()
@ -132,6 +128,12 @@ colvarproxy_lammps::~colvarproxy_lammps()
}
}
// re-initialize data where needed
void colvarproxy_lammps::setup()
{
colvars->setup();
}
// trigger colvars computation
double colvarproxy_lammps::compute()
{
@ -338,7 +340,7 @@ void colvarproxy_lammps::backup_file(char const *filename)
int colvarproxy_lammps::init_lammps_atom(const int &aid, cvm::atom *atom)
{
atom->id = aid;
atom->mass = _lmp->atom->mass[_typemap[aid]];
atom->mass = 0.0;
for (size_t i = 0; i < colvars_atoms.size(); i++) {
if (colvars_atoms[i] == aid) {
@ -353,7 +355,7 @@ int colvarproxy_lammps::init_lammps_atom(const int &aid, cvm::atom *atom)
colvars_atoms.push_back(aid);
struct commdata c;
c.tag = aid;
c.type = _typemap[aid];
c.type = 0;
c.x = c.y = c.z = 0.0;
positions.push_back(c);
total_forces.push_back(c);
@ -422,13 +424,13 @@ cvm::atom::~atom()
}
}
void cvm::atom::read_position()
{
colvarproxy_lammps const * const cp = (colvarproxy_lammps *) cvm::proxy;
this->pos.x = cp->positions[this->index].x;
this->pos.y = cp->positions[this->index].y;
this->pos.z = cp->positions[this->index].z;
this->mass = cp->positions[this->index].m;
}
void cvm::atom::read_velocity()

View File

@ -17,7 +17,7 @@
/* struct for packed data communication of coordinates and forces. */
struct commdata {
int tag,type;
double x,y,z;
double x,y,z,m;
};
inline std::ostream & operator<< (std::ostream &out, const commdata &cd)
@ -38,9 +38,6 @@ class colvarproxy_lammps : public colvarproxy {
class LAMMPS_NS::LAMMPS *_lmp;
class LAMMPS_NS::RanPark *_random;
// pointers to LAMMPS provided storage
const int *_typemap;
// state of LAMMPS properties
double t_target;
double bias_energy;
@ -66,7 +63,8 @@ class colvarproxy_lammps : public colvarproxy {
colvarproxy_lammps (LAMMPS_NS::LAMMPS *lmp, const char *,
const char *, const int, const double);
virtual ~colvarproxy_lammps();
void init(const char*, const int *);
void init(const char*);
void setup();
// disable default and copy constructor
private:

View File

@ -293,9 +293,6 @@ FixColvars::FixColvars(LAMMPS *lmp, int narg, char **arg) :
if (narg < 4)
error->all(FLERR,"Illegal fix colvars command: too few arguments");
if (atom->rmass_flag)
error->all(FLERR,"Cannot use fix colvars for atoms with rmass attribute");
if (instances > 0)
error->all(FLERR,"Only one colvars fix can be active at a time");
++instances;
@ -420,7 +417,6 @@ void FixColvars::init()
void FixColvars::setup(int vflag)
{
int *typemap,*type_buf;
const int * const tag = atom->tag;
const int * const type = atom->type;
int i,nme,tmp,ndata,nlocal_max,tag_max,max;
@ -433,6 +429,8 @@ void FixColvars::setup(int vflag)
if (init_flag == 0) {
init_flag = 1;
#if 0
// collect a list of atom type by atom id for the entire system.
// the colvar module requires this information to set masses. :-(
@ -478,6 +476,7 @@ void FixColvars::setup(int vflag)
MPI_Rsend(type_buf, nme, MPI_INT, 0, 0, world);
}
delete type_buf;
#endif
// now create and initialize the colvars proxy
@ -505,13 +504,11 @@ void FixColvars::setup(int vflag)
}
proxy = new colvarproxy_lammps(lmp,inp_name,out_name,rng_seed,t_target);
proxy->init(conf_file,typemap);
proxy->init(conf_file);
coords = proxy->get_coords();
forces = proxy->get_forces();
oforce = proxy->get_oforce();
num_coords = coords->size();
delete typemap;
}
// send the list of all colvar atom IDs to all nodes.
@ -585,6 +582,11 @@ void FixColvars::setup(int vflag)
cd[i].y = x[k][1];
cd[i].z = x[k][2];
}
if (atom->rmass_flag) {
cd[i].m = atom->rmass[k];
} else {
cd[i].m = atom->mass[type[k]];
}
}
}
@ -606,6 +608,7 @@ void FixColvars::setup(int vflag)
cd[j].x = comm_buf[k].x;
cd[j].y = comm_buf[k].y;
cd[j].z = comm_buf[k].z;
cd[j].m = comm_buf[k].m;
of[j].x = of[j].y = of[j].z = 0.0;
}
}
@ -636,6 +639,12 @@ void FixColvars::setup(int vflag)
comm_buf[nme].z = x[k][2];
}
if (atom->rmass_flag) {
comm_buf[nme].m = atom->rmass[k];
} else {
comm_buf[nme].m = atom->mass[type[k]];
}
++nme;
}
}
@ -644,6 +653,9 @@ void FixColvars::setup(int vflag)
MPI_Rsend(comm_buf, nme*size_one, MPI_BYTE, 0, 0, world);
}
// run pre-run setup in colvarproxy
proxy->setup();
// initialize forces
if (strstr(update->integrate_style,"verlet") || (update->whichflag == 2))
post_force(vflag);