1 // clang-format off
2 /* ----------------------------------------------------------------------
3    LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
4    https://www.lammps.org/
5    Steve Plimpton, sjplimp@sandia.gov, Sandia National Laboratories
6 
7    Copyright (2003) Sandia Corporation.  Under the terms of Contract
8    DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9    certain rights in this software.  This software is distributed under
10    the GNU General Public License.
11 
12    See the README file in the top-level LAMMPS directory.
13 ------------------------------------------------------------------------- */
14 
15 /* ----------------------------------------------------------------------
16    Contributing author: Pieter in 't Veld (SNL)
17 ------------------------------------------------------------------------- */
18 
19 #include "fix_nvt_sllod.h"
20 
21 #include "atom.h"
22 #include "compute.h"
23 #include "domain.h"
24 #include "error.h"
25 #include "fix_deform.h"
26 #include "group.h"
27 #include "math_extra.h"
28 #include "modify.h"
29 
30 #include <cstring>
31 
32 using namespace LAMMPS_NS;
33 using namespace FixConst;
34 
35 /* ---------------------------------------------------------------------- */
36 
FixNVTSllod(LAMMPS * lmp,int narg,char ** arg)37 FixNVTSllod::FixNVTSllod(LAMMPS *lmp, int narg, char **arg) :
38   FixNH(lmp, narg, arg)
39 {
40   if (!tstat_flag)
41     error->all(FLERR,"Temperature control must be used with fix nvt/sllod");
42   if (pstat_flag)
43     error->all(FLERR,"Pressure control can not be used with fix nvt/sllod");
44 
45   // default values
46 
47   if (mtchain_default_flag) mtchain = 1;
48 
49 
50   // create a new compute temp style
51   // id = fix-ID + temp
52 
53   id_temp = utils::strdup(std::string(id) + "_temp");
54   modify->add_compute(fmt::format("{} {} temp/deform",
55                                   id_temp,group->names[igroup]));
56   tcomputeflag = 1;
57 }
58 
59 /* ---------------------------------------------------------------------- */
60 
init()61 void FixNVTSllod::init()
62 {
63   FixNH::init();
64 
65   if (!temperature->tempbias)
66     error->all(FLERR,"Temperature for fix nvt/sllod does not have a bias");
67 
68   nondeformbias = 0;
69   if (strcmp(temperature->style,"temp/deform") != 0) nondeformbias = 1;
70 
71   // check fix deform remap settings
72 
73   int i;
74   for (i = 0; i < modify->nfix; i++)
75     if (strncmp(modify->fix[i]->style,"deform",6) == 0) {
76       if (((FixDeform *) modify->fix[i])->remapflag != Domain::V_REMAP)
77         error->all(FLERR,"Using fix nvt/sllod with inconsistent fix deform "
78                    "remap option");
79       break;
80     }
81   if (i == modify->nfix)
82     error->all(FLERR,"Using fix nvt/sllod with no fix deform defined");
83 }
84 
85 /* ----------------------------------------------------------------------
86    perform half-step scaling of velocities
87 -----------------------------------------------------------------------*/
88 
nh_v_temp()89 void FixNVTSllod::nh_v_temp()
90 {
91   // remove and restore bias = streaming velocity = Hrate*lamda + Hratelo
92   // thermostat thermal velocity only
93   // vdelu = SLLOD correction = Hrate*Hinv*vthermal
94   // for non temp/deform BIAS:
95   //   calculate temperature since some computes require temp
96   //   computed on current nlocal atoms to remove bias
97 
98   if (nondeformbias) temperature->compute_scalar();
99 
100   double **v = atom->v;
101   int *mask = atom->mask;
102   int nlocal = atom->nlocal;
103   if (igroup == atom->firstgroup) nlocal = atom->nfirst;
104 
105   double h_two[6],vdelu[3];
106   MathExtra::multiply_shape_shape(domain->h_rate,domain->h_inv,h_two);
107 
108   for (int i = 0; i < nlocal; i++) {
109     if (mask[i] & groupbit) {
110       vdelu[0] = h_two[0]*v[i][0] + h_two[5]*v[i][1] + h_two[4]*v[i][2];
111       vdelu[1] = h_two[1]*v[i][1] + h_two[3]*v[i][2];
112       vdelu[2] = h_two[2]*v[i][2];
113       temperature->remove_bias(i,v[i]);
114       v[i][0] = v[i][0]*factor_eta - dthalf*vdelu[0];
115       v[i][1] = v[i][1]*factor_eta - dthalf*vdelu[1];
116       v[i][2] = v[i][2]*factor_eta - dthalf*vdelu[2];
117       temperature->restore_bias(i,v[i]);
118     }
119   }
120 }
121