1 // clang-format off
2 /* ----------------------------------------------------------------------
3 LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
4 https://www.lammps.org/, Sandia National Laboratories
5 Steve Plimpton, sjplimp@sandia.gov
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: Georgios G. Vogiatzis (CoMSE, NTU Athens),
17 gvog@chemeng.ntua.gr
18 ------------------------------------------------------------------------- */
19
20 #include "improper_cossq.h"
21
22 #include "atom.h"
23 #include "comm.h"
24 #include "error.h"
25 #include "force.h"
26 #include "math_const.h"
27 #include "memory.h"
28 #include "neighbor.h"
29
30 #include <cmath>
31
32 using namespace LAMMPS_NS;
33 using namespace MathConst;
34
35 #define TOLERANCE 0.05
36 #define SMALL 0.001
37
38 /* ---------------------------------------------------------------------- */
39
ImproperCossq(LAMMPS * lmp)40 ImproperCossq::ImproperCossq(LAMMPS *lmp) : Improper(lmp) {}
41
42 /* ---------------------------------------------------------------------- */
43
~ImproperCossq()44 ImproperCossq::~ImproperCossq()
45 {
46 if (allocated) {
47 memory->destroy(setflag);
48 memory->destroy(k);
49 memory->destroy(chi);
50 }
51 }
52
53 /* ---------------------------------------------------------------------- */
54
compute(int eflag,int vflag)55 void ImproperCossq::compute(int eflag, int vflag)
56 {
57 int i1,i2,i3,i4,n,type;
58 double vb1x,vb1y,vb1z,vb2x,vb2y,vb2z,vb3x,vb3y,vb3z ;
59 double eimproper,f1[3],f2[3],f3[3],f4[3];
60 double rjisq, rji, rlksq, rlk, cosphi, angfac;
61 double cjiji, clkji, clklk, cfact1, cfact2, cfact3;
62
63 eimproper = 0.0;
64 ev_init(eflag,vflag);
65
66 double **x = atom->x;
67 double **f = atom->f;
68 int **improperlist = neighbor->improperlist;
69 int nimproperlist = neighbor->nimproperlist;
70 int nlocal = atom->nlocal;
71 int newton_bond = force->newton_bond;
72
73 for (n = 0; n < nimproperlist; n++) {
74 /* Ask the improper list for the atom types. */
75 i1 = improperlist[n][0];
76 i2 = improperlist[n][1];
77 i3 = improperlist[n][2];
78 i4 = improperlist[n][3];
79 type = improperlist[n][4];
80
81 /* separation vector between i1 and i2, (i2-i1) */
82 vb1x = x[i2][0] - x[i1][0];
83 vb1y = x[i2][1] - x[i1][1];
84 vb1z = x[i2][2] - x[i1][2];
85 rjisq = vb1x*vb1x + vb1y*vb1y + vb1z*vb1z ;
86 rji = sqrt(rjisq);
87
88 /* separation vector between i2 and i3 (i3-i2) */
89 vb2x = x[i3][0] - x[i2][0];
90 vb2y = x[i3][1] - x[i2][1];
91 vb2z = x[i3][2] - x[i2][2];
92
93 /* separation vector between i3 and i4, (i4-i3) */
94 vb3x = x[i4][0] - x[i3][0];
95 vb3y = x[i4][1] - x[i3][1];
96 vb3z = x[i4][2] - x[i3][2];
97 rlksq = vb3x*vb3x + vb3y*vb3y + vb3z*vb3z ;
98 rlk = sqrt(rlksq);
99
100 cosphi = (vb3x*vb1x + vb3y*vb1y + vb3z*vb1z)/(rji * rlk);
101
102 /* Check that cos(phi) is in the correct limits. */
103 if (cosphi > 1.0 + TOLERANCE || cosphi < (-1.0 - TOLERANCE))
104 problem(FLERR, i1, i2, i3, i4);
105
106 /* Apply corrections to round-off errors. */
107 if (cosphi > 1.0) cosphi -= SMALL;
108 if (cosphi < -1.0) cosphi += SMALL;
109
110 /* Calculate the angle: */
111 double torangle = acos(cosphi);
112 cosphi = cos(torangle - chi[type]);
113
114 if (eflag) eimproper = 0.5 * k[type] * cosphi * cosphi;
115
116 /*
117 printf("The tags: %d-%d-%d-%d, of type %d .\n",atom->tag[i1],atom->tag[i2],atom->tag[i3],atom->tag[i4],type);
118 printf("The ji vector: %f, %f, %f.\nThe lk vector: %f, %f, %f.\n", vb1x,vb1y,vb1z,vb3x,vb3y,vb3z);
119 printf("The cosine of the angle: %-1.16e.\n", cosphi);
120 printf("The energy of the improper: %-1.16e with prefactor %-1.16e.\n", eimproper, 0.5*k[type]);
121 */
122
123 /* Work out forces. */
124 angfac = - k[type] * cosphi;
125
126 cjiji = rjisq;
127 clklk = rlksq;
128 /*CLKJI = RXLK * RXJI + RYLK * RYJI + RZLK * RZJI */
129 clkji = vb3x*vb1x + vb3y*vb1y + vb3z*vb1z;
130
131 /*CFACT1 = CLKLK * CJIJI
132 CFACT1 = SQRT(CFACT1)
133 CFACT1 = ANGFAC / CFACT1*/
134 cfact1 = angfac / sqrt(clklk * cjiji);
135 /*CFACT2 = CLKJI / CLKLK*/
136 cfact2 = clkji / clklk;
137 /*CFACT3 = CLKJI / CJIJI*/
138 cfact3 = clkji / cjiji;
139
140 /*FIX = -RXLK + CFACT3 * RXJI
141 FIY = -RYLK + CFACT3 * RYJI
142 FIZ = -RZLK + CFACT3 * RZJI*/
143 f1[0] = - vb3x + cfact3 * vb1x;
144 f1[1] = - vb3y + cfact3 * vb1y;
145 f1[2] = - vb3z + cfact3 * vb1z;
146
147 /*FJX = -FIX
148 FJY = -FIY
149 FJZ = -FIZ*/
150 f2[0] = - f1[0];
151 f2[1] = - f1[1];
152 f2[2] = - f1[2];
153
154 /*FKX = CFACT2 * RXLK - RXJI
155 FKY = CFACT2 * RYLK - RYJI
156 FKZ = CFACT2 * RZLK - RZJI*/
157 f3[0] = cfact2 * vb3x - vb1x;
158 f3[1] = cfact2 * vb3y - vb1y;
159 f3[2] = cfact2 * vb3z - vb1z;
160
161 /*FLX = -FKX
162 FLY = -FKY
163 FLZ = -FKZ*/
164 f4[0] = - f3[0];
165 f4[1] = - f3[1];
166 f4[2] = - f3[2];
167
168 /*FIX = FIX * CFACT1
169 FIY = FIY * CFACT1
170 FIZ = FIZ * CFACT1*/
171 f1[0] *= cfact1;
172 f1[1] *= cfact1;
173 f1[2] *= cfact1;
174
175 /*FJX = FJX * CFACT1
176 FJY = FJY * CFACT1
177 FJZ = FJZ * CFACT1*/
178 f2[0] *= cfact1;
179 f2[1] *= cfact1;
180 f2[2] *= cfact1;
181
182 /*FKX = FKX * CFACT1
183 FKY = FKY * CFACT1
184 FKZ = FKZ * CFACT1*/
185 f3[0] *= cfact1;
186 f3[1] *= cfact1;
187 f3[2] *= cfact1;
188
189 /*FLX = FLX * CFACT1
190 FLY = FLY * CFACT1
191 FLZ = FLZ * CFACT1*/
192 f4[0] *= cfact1;
193 f4[1] *= cfact1;
194 f4[2] *= cfact1;
195
196 /* Apply force to each of 4 atoms */
197 if (newton_bond || i1 < nlocal) {
198 f[i1][0] += f1[0];
199 f[i1][1] += f1[1];
200 f[i1][2] += f1[2];
201 }
202
203 if (newton_bond || i2 < nlocal) {
204 f[i2][0] += f2[0];
205 f[i2][1] += f2[1];
206 f[i2][2] += f2[2];
207 }
208
209 if (newton_bond || i3 < nlocal) {
210 f[i3][0] += f3[0];
211 f[i3][1] += f3[1];
212 f[i3][2] += f3[2];
213 }
214
215 if (newton_bond || i4 < nlocal) {
216 f[i4][0] += f4[0];
217 f[i4][1] += f4[1];
218 f[i4][2] += f4[2];
219 }
220
221 if (evflag)
222 ev_tally(i1,i2,i3,i4,nlocal,newton_bond,eimproper,f1,f3,f4,
223 -vb1x,-vb1y,-vb1z,vb2x,vb2y,vb2z,vb3x,vb3y,vb3z);
224 }
225 }
226
227 /* ---------------------------------------------------------------------- */
228
allocate()229 void ImproperCossq::allocate()
230 {
231 allocated = 1;
232 int n = atom->nimpropertypes;
233
234 memory->create(k,n+1,"improper:k");
235 memory->create(chi,n+1,"improper:chi");
236
237 memory->create(setflag,n+1,"improper:setflag");
238 for (int i = 1; i <= n; i++) setflag[i] = 0;
239 }
240
241 /* ----------------------------------------------------------------------
242 set coeffs for one type
243 ------------------------------------------------------------------------- */
244
coeff(int narg,char ** arg)245 void ImproperCossq::coeff(int narg, char **arg)
246 {
247 /* Check whether there exist sufficient number of arguments.
248 0: type of improper to be applied to
249 1: energetic constant
250 2: equilibrium angle in degrees */
251 if (narg != 3) error->all(FLERR,"Incorrect args for cossq improper coefficients");
252 if (!allocated) allocate();
253
254 int ilo,ihi;
255 utils::bounds(FLERR,arg[0],1,atom->nimpropertypes,ilo,ihi,error);
256
257 double k_one = utils::numeric(FLERR,arg[1],false,lmp);
258 double chi_one = utils::numeric(FLERR,arg[2],false,lmp);
259
260 int count = 0;
261 for (int i = ilo; i <= ihi; i++) {
262 k[i] = k_one;
263 chi[i] = ((chi_one * MY_PI)/180.0);
264 setflag[i] = 1;
265 count++;
266 }
267
268 if (count == 0) error->all(FLERR,"Incorrect args for improper coefficients");
269 }
270
271 /* ----------------------------------------------------------------------
272 proc 0 writes out coeffs to restart file
273 ------------------------------------------------------------------------- */
write_restart(FILE * fp)274 void ImproperCossq::write_restart(FILE *fp)
275 {
276 fwrite(&k[1],sizeof(double),atom->nimpropertypes,fp);
277 fwrite(&chi[1],sizeof(double),atom->nimpropertypes,fp);
278 }
279
280 /* ----------------------------------------------------------------------
281 proc 0 reads coeffs from restart file, bcasts them
282 ------------------------------------------------------------------------- */
read_restart(FILE * fp)283 void ImproperCossq::read_restart(FILE *fp)
284 {
285 allocate();
286
287 if (comm->me == 0) {
288 utils::sfread(FLERR,&k[1],sizeof(double),atom->nimpropertypes,fp,nullptr,error);
289 utils::sfread(FLERR,&chi[1],sizeof(double),atom->nimpropertypes,fp,nullptr,error);
290 }
291 MPI_Bcast(&k[1],atom->nimpropertypes,MPI_DOUBLE,0,world);
292 MPI_Bcast(&chi[1],atom->nimpropertypes,MPI_DOUBLE,0,world);
293
294 for (int i = 1; i <= atom->nimpropertypes; i++) setflag[i] = 1;
295 }
296
297 /* ----------------------------------------------------------------------
298 proc 0 writes to data file
299 ------------------------------------------------------------------------- */
300
write_data(FILE * fp)301 void ImproperCossq::write_data(FILE *fp)
302 {
303 for (int i = 1; i <= atom->nimpropertypes; i++)
304 fprintf(fp,"%d %g %g\n",i,k[i],chi[i]/MY_PI*180.0);
305 }
306