Added FixMSST in package shock

git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@4275 f3b2605a-c512-4ea7-a41b-209d697bcdaa
This commit is contained in:
athomps 2010-06-16 21:05:27 +00:00
parent d9ce487090
commit 73b32571d9
5 changed files with 1273 additions and 0 deletions

19
src/SHOCK/Install.sh Normal file
View File

@ -0,0 +1,19 @@
# Install/unInstall package files in LAMMPS
if (test $1 = 1) then
cp fix_msst.cpp ..
cp compute_vsum.cpp ..
cp fix_msst.h ..
cp compute_vsum.h ..
elif (test $1 = 0) then
rm ../fix_msst.cpp
rm ../compute_vsum.cpp
rm ../fix_msst.h
rm ../compute_vsum.h
fi

View File

@ -0,0 +1,77 @@
/* ----------------------------------------------------------------------
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 "mpi.h"
#include "compute_vsum.h"
#include "atom.h"
#include "force.h"
#include "modify.h"
#include "fix.h"
#include "group.h"
#include "error.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
ComputeVsum::ComputeVsum(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg != 3) error->all("Illegal compute temp command");
scalar_flag = vector_flag = 1;
size_vector = 6;
extscalar = 0;
tempflag = 1;
vector = new double[6];
}
/* ---------------------------------------------------------------------- */
ComputeVsum::~ComputeVsum()
{
delete [] vector;
}
/* ---------------------------------------------------------------------- */
void ComputeVsum::init()
{
;
}
/* ---------------------------------------------------------------------- */
/* ---------------------------------------------------------------------- */
double ComputeVsum::compute_scalar()
{
double **v = atom->v;
int *mask = atom->mask;
int nlocal = atom->nlocal;
double t = 0.0;
for (int i = 0; i < nlocal; i++) {
if (mask[i] & groupbit) {
t += (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]) ;
}
}
MPI_Allreduce(&t,&scalar,1,MPI_DOUBLE,MPI_SUM,world);
return scalar;
}
/* ---------------------------------------------------------------------- */

42
src/SHOCK/compute_vsum.h Normal file
View File

@ -0,0 +1,42 @@
/* ----------------------------------------------------------------------
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.
------------------------------------------------------------------------- */
// Compute the sum of the squared velocities, without usual mass terms
// for the kinetic energy. This is used for the MSST (Larry Fried).
#ifdef COMPUTE_CLASS
ComputeStyle(vsum,ComputeVsum)
#else
#ifndef COMPUTE_VSUM_H
#define COMPUTE_VSUM_H
#include "compute.h"
namespace LAMMPS_NS {
class ComputeVsum : public Compute {
public:
ComputeVsum(class LAMMPS *, int, char **);
~ComputeVsum();
void init();
double compute_scalar();
};
}
#endif
#endif

1027
src/SHOCK/fix_msst.cpp Normal file

File diff suppressed because it is too large Load Diff

108
src/SHOCK/fix_msst.h Normal file
View File

@ -0,0 +1,108 @@
/* ----------------------------------------------------------------------
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.
------------------------------------------------------------------------- */
/* Implementation of the Multi-Scale Shock Method.
See Reed, Fried, Joannopoulos, Phys. Rev. Lett., 90, 235503(2003).
Implementation by Laurence Fried, LLNL, 4/2007.
*/
#ifdef FIX_CLASS
FixStyle(msst,FixMSST)
#else
#ifndef FIX_MSST_H
#define FIX_MSST_H
#include "fix.h"
namespace LAMMPS_NS {
class FixMSST : public Fix {
public:
FixMSST(class LAMMPS *, int, char **);
~FixMSST();
int setmask();
void init();
void setup(int);
void initial_integrate(int);
void final_integrate();
void initial_integrate_respa(int,int,int);
void final_integrate_respa(int);
double compute_scalar();
double compute_vector(int);
void write_restart(FILE *);
void restart(char *);
int modify_param(int, char **);
private:
double dtv,dtf,dthalf; // Full and half step sizes.
double boltz,nktv2p, mvv2e; // Boltzmann factor and unit conversions.
double total_mass; // Mass of the computational cell.
double omega[3]; // Time derivative of the volume.
double p_current[3],dilation[3];
double qmass; // Effective cell mass.
double mu; // Effective cell viscosity.
double tscale; // Converts thermal energy to compressive
// strain ke at simulation start
double velocity_sum; // Sum of the velocities squared.
double **old_velocity; // Saved velocities.
int kspace_flag; // 1 if KSpace invoked, 0 if not
int nrigid; // number of rigid fixes
int *rfix; // indices of rigid fixes
char *id_temp,*id_press; // Strings with identifiers of
char *id_vsum, *id_pe; // created computes.
class Compute *temperature; // Computes created to evaluate
class Compute *pressure; // thermodynamic quantities.
class Compute *vsum;
class Compute *pe;
int tflag,pflag,vsflag,peflag; // Flags to keep track of computes that
// were created.
// shock initial conditions.
double e0; // Initial energy
double v0; // Initial volume
double p0; // Initial pressure
double velocity; // Velocity of the shock.
double lagrangian_position; // Lagrangian location of computational cell
int direction; // Direction of shock
int p0_set; // Is pressure set.
int v0_set; // Is volume set.
int e0_set; // Is energy set.
int atoms_allocated; // The number of allocated atoms in old_velocity.
// functions
void couple();
void remap(int);
void check_alloc(int n);
double compute_etotal();
double compute_vol();
double compute_hugoniot();
double compute_rayleigh();
double compute_lagrangian_speed();
double compute_lagrangian_position();
};
}
#endif
#endif