/* ---------------------------------------------------------------------- LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator http://lammps.sandia.gov, Sandia National Laboratories Steve Plimpton, sjplimp@sandia.gov Copyright (2003) Sandia Corporation. Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains certain rights in this software. This software is distributed under the GNU General Public License. See the README file in the top-level LAMMPS directory. ------------------------------------------------------------------------- */ #include #include #include "min_quickmin.h" #include "universe.h" #include "atom.h" #include "force.h" #include "update.h" #include "output.h" #include "timer.h" #include "error.h" using namespace LAMMPS_NS; // EPS_ENERGY = minimum normalization for energy tolerance #define EPS_ENERGY 1.0e-8 #define DELAYSTEP 5 /* ---------------------------------------------------------------------- */ MinQuickMin::MinQuickMin(LAMMPS *lmp) : Min(lmp) {} /* ---------------------------------------------------------------------- */ void MinQuickMin::init() { Min::init(); dt = update->dt; last_negative = update->ntimestep; } /* ---------------------------------------------------------------------- */ void MinQuickMin::setup_style() { double **v = atom->v; int nlocal = atom->nlocal; for (int i = 0; i < nlocal; i++) v[i][0] = v[i][1] = v[i][2] = 0.0; } /* ---------------------------------------------------------------------- set current vector lengths and pointers called after atoms have migrated ------------------------------------------------------------------------- */ void MinQuickMin::reset_vectors() { // atomic dof nvec = 3 * atom->nlocal; if (nvec) xvec = atom->x[0]; if (nvec) fvec = atom->f[0]; } /* ---------------------------------------------------------------------- minimization via QuickMin damped dynamics ------------------------------------------------------------------------- */ int MinQuickMin::iterate(int maxiter) { bigint ntimestep; double vmax,vdotf,vdotfall,fdotf,fdotfall,scale; double dtvone,dtv,dtf,dtfm; int flag,flagall; alpha_final = 0.0; for (int iter = 0; iter < maxiter; iter++) { if (timer->check_timeout(niter)) return TIMEOUT; ntimestep = ++update->ntimestep; niter++; // zero velocity if anti-parallel to force // else project velocity in direction of force double **v = atom->v; double **f = atom->f; int nlocal = atom->nlocal; vdotf = 0.0; for (int i = 0; i < nlocal; i++) vdotf += v[i][0]*f[i][0] + v[i][1]*f[i][1] + v[i][2]*f[i][2]; MPI_Allreduce(&vdotf,&vdotfall,1,MPI_DOUBLE,MPI_SUM,world); // sum vdotf over replicas, if necessary // this communicator would be invalid for multiprocess replicas if (update->multireplica == 1) { vdotf = vdotfall; MPI_Allreduce(&vdotf,&vdotfall,1,MPI_DOUBLE,MPI_SUM,universe->uworld); } if (vdotfall < 0.0) { last_negative = ntimestep; for (int i = 0; i < nlocal; i++) v[i][0] = v[i][1] = v[i][2] = 0.0; } else { fdotf = 0.0; for (int i = 0; i < nlocal; i++) fdotf += f[i][0]*f[i][0] + f[i][1]*f[i][1] + f[i][2]*f[i][2]; MPI_Allreduce(&fdotf,&fdotfall,1,MPI_DOUBLE,MPI_SUM,world); // sum fdotf over replicas, if necessary // this communicator would be invalid for multiprocess replicas if (update->multireplica == 1) { fdotf = fdotfall; MPI_Allreduce(&fdotf,&fdotfall,1,MPI_DOUBLE,MPI_SUM,universe->uworld); } if (fdotfall == 0.0) scale = 0.0; else scale = vdotfall/fdotfall; for (int i = 0; i < nlocal; i++) { v[i][0] = scale * f[i][0]; v[i][1] = scale * f[i][1]; v[i][2] = scale * f[i][2]; } } // limit timestep so no particle moves further than dmax double *rmass = atom->rmass; double *mass = atom->mass; int *type = atom->type; dtvone = dt; for (int i = 0; i < nlocal; i++) { vmax = MAX(fabs(v[i][0]),fabs(v[i][1])); vmax = MAX(vmax,fabs(v[i][2])); if (dtvone*vmax > dmax) dtvone = dmax/vmax; } MPI_Allreduce(&dtvone,&dtv,1,MPI_DOUBLE,MPI_MIN,world); // min dtv over replicas, if necessary // this communicator would be invalid for multiprocess replicas if (update->multireplica == 1) { dtvone = dtv; MPI_Allreduce(&dtvone,&dtv,1,MPI_DOUBLE,MPI_MIN,universe->uworld); } dtf = dtv * force->ftm2v; // Euler integration step double **x = atom->x; if (rmass) { for (int i = 0; i < nlocal; i++) { dtfm = dtf / rmass[i]; x[i][0] += dtv * v[i][0]; x[i][1] += dtv * v[i][1]; x[i][2] += dtv * v[i][2]; v[i][0] += dtfm * f[i][0]; v[i][1] += dtfm * f[i][1]; v[i][2] += dtfm * f[i][2]; } } else { for (int i = 0; i < nlocal; i++) { dtfm = dtf / mass[type[i]]; x[i][0] += dtv * v[i][0]; x[i][1] += dtv * v[i][1]; x[i][2] += dtv * v[i][2]; v[i][0] += dtfm * f[i][0]; v[i][1] += dtfm * f[i][1]; v[i][2] += dtfm * f[i][2]; } } eprevious = ecurrent; ecurrent = energy_force(0); neval++; // energy tolerance criterion // only check after DELAYSTEP elapsed since velocties reset to 0 // sync across replicas if running multi-replica minimization if (update->etol > 0.0 && ntimestep-last_negative > DELAYSTEP) { if (update->multireplica == 0) { if (fabs(ecurrent-eprevious) < update->etol * 0.5*(fabs(ecurrent) + fabs(eprevious) + EPS_ENERGY)) return ETOL; } else { if (fabs(ecurrent-eprevious) < update->etol * 0.5*(fabs(ecurrent) + fabs(eprevious) + EPS_ENERGY)) flag = 0; else flag = 1; MPI_Allreduce(&flag,&flagall,1,MPI_INT,MPI_SUM,universe->uworld); if (flagall == 0) return ETOL; } } // force tolerance criterion // sync across replicas if running multi-replica minimization if (update->ftol > 0.0) { fdotf = fnorm_sqr(); if (update->multireplica == 0) { if (fdotf < update->ftol*update->ftol) return FTOL; } else { if (fdotf < update->ftol*update->ftol) flag = 0; else flag = 1; MPI_Allreduce(&flag,&flagall,1,MPI_INT,MPI_SUM,universe->uworld); if (flagall == 0) return FTOL; } } // output for thermo, dump, restart files if (output->next == ntimestep) { timer->stamp(); output->write(ntimestep); timer->stamp(Timer::OUTPUT); } } return MAXITER; }