forked from lammps/lammps
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathfix_thermal_conductivity.cpp
328 lines (276 loc) · 10.6 KB
/
fix_thermal_conductivity.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, [email protected]
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.
------------------------------------------------------------------------- */
/* ----------------------------------------------------------------------
Contributing author: Craig Tenney (UND) added support
for swapping atoms of different masses
------------------------------------------------------------------------- */
#include "fix_thermal_conductivity.h"
#include <cstring>
#include "atom.h"
#include "force.h"
#include "domain.h"
#include "modify.h"
#include "error.h"
using namespace LAMMPS_NS;
using namespace FixConst;
#define BIG 1.0e10
/* ---------------------------------------------------------------------- */
FixThermalConductivity::FixThermalConductivity(LAMMPS *lmp,
int narg, char **arg) :
Fix(lmp, narg, arg),
index_lo(nullptr), index_hi(nullptr), ke_lo(nullptr), ke_hi(nullptr)
{
if (narg < 6) error->all(FLERR,"Illegal fix thermal/conductivity command");
MPI_Comm_rank(world,&me);
nevery = utils::inumeric(FLERR,arg[3],false,lmp);
if (nevery <= 0) error->all(FLERR,"Illegal fix thermal/conductivity command");
scalar_flag = 1;
global_freq = nevery;
extscalar = 0;
if (strcmp(arg[4],"x") == 0) edim = 0;
else if (strcmp(arg[4],"y") == 0) edim = 1;
else if (strcmp(arg[4],"z") == 0) edim = 2;
else error->all(FLERR,"Illegal fix thermal/conductivity command");
nbin = utils::inumeric(FLERR,arg[5],false,lmp);
if (nbin % 2 || nbin <= 2)
error->all(FLERR,"Illegal fix thermal/conductivity command");
// optional keywords
nswap = 1;
int iarg = 6;
while (iarg < narg) {
if (strcmp(arg[iarg],"swap") == 0) {
if (iarg+2 > narg)
error->all(FLERR,"Illegal fix thermal/conductivity command");
nswap = utils::inumeric(FLERR,arg[iarg+1],false,lmp);
if (nswap <= 0)
error->all(FLERR,
"Fix thermal/conductivity swap value must be positive");
iarg += 2;
} else error->all(FLERR,"Illegal fix thermal/conductivity command");
}
// initialize array sizes to nswap+1 so have space to shift values down
index_lo = new int[nswap+1];
index_hi = new int[nswap+1];
ke_lo = new double[nswap+1];
ke_hi = new double[nswap+1];
e_exchange = 0.0;
}
/* ---------------------------------------------------------------------- */
FixThermalConductivity::~FixThermalConductivity()
{
delete [] index_lo;
delete [] index_hi;
delete [] ke_lo;
delete [] ke_hi;
}
/* ---------------------------------------------------------------------- */
int FixThermalConductivity::setmask()
{
int mask = 0;
mask |= END_OF_STEP;
return mask;
}
/* ---------------------------------------------------------------------- */
void FixThermalConductivity::init()
{
// warn if any fix ave/spatial comes after this fix
// can cause glitch in averaging since ave will happen after swap
int foundme = 0;
for (int i = 0; i < modify->nfix; i++) {
if (modify->fix[i] == this) foundme = 1;
if (foundme && strcmp(modify->fix[i]->style,"ave/spatial") == 0 && me == 0)
error->warning(FLERR,
"Fix thermal/conductivity comes before fix ave/spatial");
}
// set bounds of 2 slabs in edim
// only necessary for static box, else re-computed in end_of_step()
// lo bin is always bottom bin
// hi bin is just above half height
if (domain->box_change == 0) {
prd = domain->prd[edim];
boxlo = domain->boxlo[edim];
boxhi = domain->boxhi[edim];
double binsize = (boxhi-boxlo) / nbin;
slablo_lo = boxlo;
slablo_hi = boxlo + binsize;
slabhi_lo = boxlo + (nbin/2)*binsize;
slabhi_hi = boxlo + (nbin/2+1)*binsize;
}
periodicity = domain->periodicity[edim];
}
/* ---------------------------------------------------------------------- */
void FixThermalConductivity::end_of_step()
{
int i,j,m,insert;
double coord,ke;
struct {
double value;
int proc;
} mine[2],all[2];
// if box changes, recompute bounds of 2 slabs in edim
if (domain->box_change) {
prd = domain->prd[edim];
boxlo = domain->boxlo[edim];
boxhi = domain->boxhi[edim];
double binsize = (boxhi-boxlo) / nbin;
slablo_lo = boxlo;
slablo_hi = boxlo + binsize;
slabhi_lo = boxlo + (nbin/2)*binsize;
slabhi_hi = boxlo + (nbin/2+1)*binsize;
}
// make 2 lists of up to nswap atoms
// hottest atoms in lo slab, coldest atoms in hi slab (really mid slab)
// lo slab list is sorted by hottest, hi slab is sorted by coldest
// map atoms back into periodic box if necessary
// insert = location in list to insert new atom
double **x = atom->x;
double **v = atom->v;
double *mass = atom->mass;
double *rmass = atom->rmass;
int *type = atom->type;
int *mask = atom->mask;
int nlocal = atom->nlocal;
nhi = nlo = 0;
for (i = 0; i < nlocal; i++)
if (mask[i] & groupbit) {
coord = x[i][edim];
if (coord < boxlo && periodicity) coord += prd;
else if (coord >= boxhi && periodicity) coord -= prd;
if (coord >= slablo_lo && coord < slablo_hi) {
ke = v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2];
if (rmass) ke *= 0.5*rmass[i];
else ke *= 0.5*mass[type[i]];
if (nlo < nswap || ke > ke_lo[nswap-1]) {
for (insert = nlo-1; insert >= 0; insert--)
if (ke < ke_lo[insert]) break;
insert++;
for (m = nlo-1; m >= insert; m--) {
ke_lo[m+1] = ke_lo[m];
index_lo[m+1] = index_lo[m];
}
ke_lo[insert] = ke;
index_lo[insert] = i;
if (nlo < nswap) nlo++;
}
}
if (coord >= slabhi_lo && coord < slabhi_hi) {
ke = v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2];
if (rmass) ke *= 0.5*rmass[i];
else ke *= 0.5*mass[type[i]];
if (nhi < nswap || ke < ke_hi[nswap-1]) {
for (insert = nhi-1; insert >= 0; insert--)
if (ke > ke_hi[insert]) break;
insert++;
for (m = nhi-1; m >= insert; m--) {
ke_hi[m+1] = ke_hi[m];
index_hi[m+1] = index_hi[m];
}
ke_hi[insert] = ke;
index_hi[insert] = i;
if (nhi < nswap) nhi++;
}
}
}
// loop over nswap pairs
// pair 2 global atoms at beginning of sorted lo/hi slab lists via Allreduce
// BIG values are for procs with no atom to contribute
// use negative of hottest KE since is doing a MINLOC
// MINLOC also communicates which procs own them
// exchange kinetic energy between the 2 particles
// if I own both particles just swap, else point2point comm of velocities
double sbuf[4],rbuf[4],vcm[3];
double eswap = 0.0;
mine[0].proc = mine[1].proc = me;
int ilo = 0;
int ihi = 0;
for (m = 0; m < nswap; m++) {
if (ilo < nlo) mine[0].value = -ke_lo[ilo];
else mine[0].value = BIG;
if (ihi < nhi) mine[1].value = ke_hi[ihi];
else mine[1].value = BIG;
MPI_Allreduce(mine,all,2,MPI_DOUBLE_INT,MPI_MINLOC,world);
if (all[0].value == BIG || all[1].value == BIG) continue;
if (me == all[0].proc && me == all[1].proc) {
i = index_lo[ilo++];
j = index_hi[ihi++];
sbuf[0] = v[j][0];
sbuf[1] = v[j][1];
sbuf[2] = v[j][2];
if (rmass) sbuf[3] = rmass[j];
else sbuf[3] = mass[type[j]];
rbuf[0] = v[i][0];
rbuf[1] = v[i][1];
rbuf[2] = v[i][2];
if (rmass) rbuf[3] = rmass[i];
else rbuf[3] = mass[type[i]];
vcm[0] = (sbuf[3]*sbuf[0] + rbuf[3]*rbuf[0]) / (sbuf[3] + rbuf[3]);
vcm[1] = (sbuf[3]*sbuf[1] + rbuf[3]*rbuf[1]) / (sbuf[3] + rbuf[3]);
vcm[2] = (sbuf[3]*sbuf[2] + rbuf[3]*rbuf[2]) / (sbuf[3] + rbuf[3]);
v[j][0] = 2.0 * vcm[0] - sbuf[0];
v[j][1] = 2.0 * vcm[1] - sbuf[1];
v[j][2] = 2.0 * vcm[2] - sbuf[2];
eswap += sbuf[3] * (vcm[0] * (vcm[0] - sbuf[0]) +
vcm[1] * (vcm[1] - sbuf[1]) +
vcm[2] * (vcm[2] - sbuf[2]));
v[i][0] = 2.0 * vcm[0] - rbuf[0];
v[i][1] = 2.0 * vcm[1] - rbuf[1];
v[i][2] = 2.0 * vcm[2] - rbuf[2];
eswap -= rbuf[3] * (vcm[0] * (vcm[0] - rbuf[0]) +
vcm[1] * (vcm[1] - rbuf[1]) +
vcm[2] * (vcm[2] - rbuf[2]));
} else if (me == all[0].proc) {
j = index_lo[ilo++];
sbuf[0] = v[j][0];
sbuf[1] = v[j][1];
sbuf[2] = v[j][2];
if (rmass) sbuf[3] = rmass[j];
else sbuf[3] = mass[type[j]];
MPI_Sendrecv(sbuf,4,MPI_DOUBLE,all[1].proc,0,
rbuf,4,MPI_DOUBLE,all[1].proc,0,world,MPI_STATUS_IGNORE);
vcm[0] = (sbuf[3]*sbuf[0] + rbuf[3]*rbuf[0]) / (sbuf[3] + rbuf[3]);
vcm[1] = (sbuf[3]*sbuf[1] + rbuf[3]*rbuf[1]) / (sbuf[3] + rbuf[3]);
vcm[2] = (sbuf[3]*sbuf[2] + rbuf[3]*rbuf[2]) / (sbuf[3] + rbuf[3]);
v[j][0] = 2.0 * vcm[0] - sbuf[0];
v[j][1] = 2.0 * vcm[1] - sbuf[1];
v[j][2] = 2.0 * vcm[2] - sbuf[2];
eswap -= sbuf[3] * (vcm[0] * (vcm[0] - sbuf[0]) +
vcm[1] * (vcm[1] - sbuf[1]) +
vcm[2] * (vcm[2] - sbuf[2]));
} else if (me == all[1].proc) {
j = index_hi[ihi++];
sbuf[0] = v[j][0];
sbuf[1] = v[j][1];
sbuf[2] = v[j][2];
if (rmass) sbuf[3] = rmass[j];
else sbuf[3] = mass[type[j]];
MPI_Sendrecv(sbuf,4,MPI_DOUBLE,all[0].proc,0,
rbuf,4,MPI_DOUBLE,all[0].proc,0,world,MPI_STATUS_IGNORE);
vcm[0] = (sbuf[3]*sbuf[0] + rbuf[3]*rbuf[0]) / (sbuf[3] + rbuf[3]);
vcm[1] = (sbuf[3]*sbuf[1] + rbuf[3]*rbuf[1]) / (sbuf[3] + rbuf[3]);
vcm[2] = (sbuf[3]*sbuf[2] + rbuf[3]*rbuf[2]) / (sbuf[3] + rbuf[3]);
v[j][0] = 2.0 * vcm[0] - sbuf[0];
v[j][1] = 2.0 * vcm[1] - sbuf[1];
v[j][2] = 2.0 * vcm[2] - sbuf[2];
eswap += sbuf[3] * (vcm[0] * (vcm[0] - sbuf[0]) +
vcm[1] * (vcm[1] - sbuf[1]) +
vcm[2] * (vcm[2] - sbuf[2]));
}
}
// tally energy exchange from all swaps
double eswap_all;
MPI_Allreduce(&eswap,&eswap_all,1,MPI_DOUBLE,MPI_SUM,world);
e_exchange += force->mvv2e * eswap_all;
}
/* ---------------------------------------------------------------------- */
double FixThermalConductivity::compute_scalar()
{
return e_exchange;
}