Files
lammps/src/min_cg.cpp
Axel Kohlmeyer 48aa45f7f0 Merge branch 'master' into lammps-icms
Resolved Conflicts (i hope):
	src/REPLICA/neb.cpp
	src/REPLICA/prd.cpp
	src/REPLICA/tad.cpp
	src/REPLICA/verlet_split.cpp
	src/USER-CUDA/atom_vec_angle_cuda.h
	src/USER-CUDA/atom_vec_atomic_cuda.h
	src/USER-CUDA/atom_vec_charge_cuda.h
	src/USER-CUDA/atom_vec_full_cuda.h
	src/USER-CUDA/compute_pe_cuda.h
	src/USER-CUDA/compute_pressure_cuda.h
	src/USER-CUDA/compute_temp_cuda.h
	src/USER-CUDA/compute_temp_partial_cuda.h
	src/USER-CUDA/cuda.cpp
	src/USER-CUDA/cuda.h
	src/USER-CUDA/cuda_data.h
	src/USER-CUDA/cuda_neigh_list.h
	src/USER-CUDA/fft3d_cuda.h
	src/USER-CUDA/fft3d_wrap_cuda.h
	src/USER-CUDA/fix_addforce_cuda.h
	src/USER-CUDA/fix_aveforce_cuda.h
	src/USER-CUDA/fix_enforce2d_cuda.h
	src/USER-CUDA/fix_freeze_cuda.h
	src/USER-CUDA/fix_gravity_cuda.h
	src/USER-CUDA/fix_nve_cuda.h
	src/USER-CUDA/fix_set_force_cuda.h
	src/USER-CUDA/fix_temp_berendsen_cuda.h
	src/USER-CUDA/fix_temp_rescale_cuda.h
	src/USER-CUDA/fix_temp_rescale_limit_cuda.h
	src/USER-CUDA/fix_viscous_cuda.h
	src/USER-CUDA/pair_born_coul_long_cuda.h
	src/USER-CUDA/pair_buck_coul_cut_cuda.h
	src/USER-CUDA/pair_buck_coul_long_cuda.h
	src/USER-CUDA/pair_buck_cuda.h
	src/USER-CUDA/pair_cg_cmm_coul_cut_cuda.cpp
	src/USER-CUDA/pair_cg_cmm_coul_cut_cuda.h
	src/USER-CUDA/pair_cg_cmm_coul_debye_cuda.cpp
	src/USER-CUDA/pair_cg_cmm_coul_debye_cuda.h
	src/USER-CUDA/pair_cg_cmm_coul_long_cuda.cpp
	src/USER-CUDA/pair_cg_cmm_coul_long_cuda.h
	src/USER-CUDA/pair_cg_cmm_cuda.cpp
	src/USER-CUDA/pair_cg_cmm_cuda.h
	src/USER-CUDA/pair_eam_cuda.h
	src/USER-CUDA/pair_gran_hooke_cuda.h
	src/USER-CUDA/pair_lj96_cut_cuda.h
	src/USER-CUDA/pair_lj_charmm_coul_charmm_cuda.h
	src/USER-CUDA/pair_lj_charmm_coul_charmm_implicit_cuda.h
	src/USER-CUDA/pair_lj_charmm_coul_long_cuda.h
	src/USER-CUDA/pair_lj_class2_coul_cut_cuda.h
	src/USER-CUDA/pair_lj_class2_coul_long_cuda.h
	src/USER-CUDA/pair_lj_class2_cuda.h
	src/USER-CUDA/pair_lj_cut_coul_cut_cuda.h
	src/USER-CUDA/pair_lj_cut_coul_debye_cuda.h
	src/USER-CUDA/pair_lj_cut_coul_long_cuda.h
	src/USER-CUDA/pair_lj_cut_cuda.h
	src/USER-CUDA/pair_lj_cut_experimental_cuda.h
	src/USER-CUDA/pair_lj_expand_cuda.h
	src/USER-CUDA/pair_lj_gromacs_coul_gromacs_cuda.h
	src/USER-CUDA/pair_lj_gromacs_cuda.h
	src/USER-CUDA/pair_lj_sdk_coul_long_cuda.cpp
	src/USER-CUDA/pair_lj_sdk_coul_long_cuda.h
	src/USER-CUDA/pair_lj_sdk_cuda.h
	src/USER-CUDA/pair_lj_smooth_cuda.h
	src/USER-CUDA/pair_morse_cuda.h
	src/USER-CUDA/pair_sw_cuda.h
	src/USER-CUDA/pair_tersoff_cuda.h
	src/USER-CUDA/pair_tersoff_zbl_cuda.h
	src/USER-CUDA/pppm_cuda.h
	src/USER-CUDA/verlet_cuda.cpp
	src/USER-CUDA/verlet_cuda.h
	src/comm.cpp
	src/finish.cpp
	src/fix_rigid_nvt.h
	src/respa.cpp
	src/run.cpp
2012-06-07 16:42:38 -04:00

185 lines
5.2 KiB
C++

/* ----------------------------------------------------------------------
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 "lmptype.h"
#include "mpi.h"
#include "math.h"
#include "string.h"
#include "min_cg.h"
#include "atom.h"
#include "update.h"
#include "output.h"
#include "timer.h"
using namespace LAMMPS_NS;
// EPS_ENERGY = minimum normalization for energy tolerance
#define EPS_ENERGY 1.0e-8
// same as in other min classes
enum{MAXITER,MAXEVAL,ETOL,FTOL,DOWNHILL,ZEROALPHA,ZEROFORCE,ZEROQUAD};
/* ---------------------------------------------------------------------- */
MinCG::MinCG(LAMMPS *lmp) : MinLineSearch(lmp) {}
/* ----------------------------------------------------------------------
minimization via conjugate gradient iterations
------------------------------------------------------------------------- */
int MinCG::iterate(int maxiter)
{
int i,m,n,fail,ntimestep;
double beta,gg,dot[2],dotall[2];
double *fatom,*gatom,*hatom;
// nlimit = max # of CG iterations before restarting
// set to ndoftotal unless too big
int nlimit = static_cast<int> (MIN(MAXSMALLINT,ndoftotal));
// initialize working vectors
for (i = 0; i < nvec; i++) h[i] = g[i] = fvec[i];
if (nextra_atom)
for (m = 0; m < nextra_atom; m++) {
fatom = fextra_atom[m];
gatom = gextra_atom[m];
hatom = hextra_atom[m];
n = extra_nlen[m];
for (i = 0; i < n; i++) hatom[i] = gatom[i] = fatom[i];
}
if (nextra_global)
for (i = 0; i < nextra_global; i++) hextra[i] = gextra[i] = fextra[i];
gg = fnorm_sqr();
for (int iter = 0; iter < maxiter; iter++) {
ntimestep = ++update->ntimestep;
niter++;
// line minimization along direction h from current atom->x
eprevious = ecurrent;
fail = (this->*linemin)(ecurrent,alpha_final);
if (fail) return fail;
// function evaluation criterion
if (neval >= update->max_eval) return MAXEVAL;
// energy tolerance criterion
if (fabs(ecurrent-eprevious) <
update->etol * 0.5*(fabs(ecurrent) + fabs(eprevious) + EPS_ENERGY))
return ETOL;
// force tolerance criterion
dot[0] = dot[1] = 0.0;
for (i = 0; i < nvec; i++) {
dot[0] += fvec[i]*fvec[i];
dot[1] += fvec[i]*g[i];
}
if (nextra_atom)
for (m = 0; m < nextra_atom; m++) {
fatom = fextra_atom[m];
gatom = gextra_atom[m];
n = extra_nlen[m];
for (i = 0; i < n; i++) {
dot[0] += fatom[i]*fatom[i];
dot[1] += fatom[i]*gatom[i];
}
}
MPI_Allreduce(dot,dotall,2,MPI_DOUBLE,MPI_SUM,world);
if (nextra_global)
for (i = 0; i < nextra_global; i++) {
dotall[0] += fextra[i]*fextra[i];
dotall[1] += fextra[i]*gextra[i];
}
if (dotall[0] < update->ftol*update->ftol) return FTOL;
// update new search direction h from new f = -Grad(x) and old g
// this is Polak-Ribieri formulation
// beta = dotall[0]/gg would be Fletcher-Reeves
// reinitialize CG every ndof iterations by setting beta = 0.0
beta = MAX(0.0,(dotall[0] - dotall[1])/gg);
if ((niter+1) % nlimit == 0) beta = 0.0;
gg = dotall[0];
for (i = 0; i < nvec; i++) {
g[i] = fvec[i];
h[i] = g[i] + beta*h[i];
}
if (nextra_atom)
for (m = 0; m < nextra_atom; m++) {
fatom = fextra_atom[m];
gatom = gextra_atom[m];
hatom = hextra_atom[m];
n = extra_nlen[m];
for (i = 0; i < n; i++) {
gatom[i] = fatom[i];
hatom[i] = gatom[i] + beta*hatom[i];
}
}
if (nextra_global)
for (i = 0; i < nextra_global; i++) {
gextra[i] = fextra[i];
hextra[i] = gextra[i] + beta*hextra[i];
}
// reinitialize CG if new search direction h is not downhill
dot[0] = 0.0;
for (i = 0; i < nvec; i++) dot[0] += g[i]*h[i];
if (nextra_atom)
for (m = 0; m < nextra_atom; m++) {
gatom = gextra_atom[m];
hatom = hextra_atom[m];
n = extra_nlen[m];
for (i = 0; i < n; i++) dot[0] += gatom[i]*hatom[i];
}
MPI_Allreduce(dot,dotall,1,MPI_DOUBLE,MPI_SUM,world);
if (nextra_global)
for (i = 0; i < nextra_global; i++)
dotall[0] += gextra[i]*hextra[i];
if (dotall[0] <= 0.0) {
for (i = 0; i < nvec; i++) h[i] = g[i];
if (nextra_atom)
for (m = 0; m < nextra_atom; m++) {
gatom = gextra_atom[m];
hatom = hextra_atom[m];
n = extra_nlen[m];
for (i = 0; i < n; i++) hatom[i] = gatom[i];
}
if (nextra_global)
for (i = 0; i < nextra_global; i++) hextra[i] = gextra[i];
}
// output for thermo, dump, restart files
if (output->next == ntimestep) {
timer->stamp();
output->write(ntimestep);
timer->stamp(Timer::OUTPUT);
}
}
return MAXITER;
}