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 authors: Frances Mackay, Santtu Ollila, Colin Denniston (UWO)
17
18 Based on fix_momentum,
19 Contributing author: Naveen Michaud-Agrawal (Johns Hopkins U)
20 ------------------------------------------------------------------------- */
21
22 #include "fix_lb_momentum.h"
23
24
25 #include <cstring>
26 #include "atom.h"
27 #include "group.h"
28 #include "error.h"
29 #include "fix_lb_fluid.h"
30 #include "modify.h"
31
32 using namespace LAMMPS_NS;
33 using namespace FixConst;
34
35
36 /* ---------------------------------------------------------------------- */
37
FixLbMomentum(LAMMPS * lmp,int narg,char ** arg)38 FixLbMomentum::FixLbMomentum(LAMMPS *lmp, int narg, char **arg) :
39 Fix(lmp, narg, arg)
40 {
41 if (narg < 4) error->all(FLERR,"Illegal fix lb/momentum command");
42 nevery = atoi(arg[3]);
43 if (nevery <= 0) error->all(FLERR,"Illegal fix lb/momentum command");
44
45 linear = 1;
46 xflag = 1;
47 yflag = 1;
48 zflag = 1;
49
50 int iarg = 4;
51 while (iarg < narg) {
52 if (strcmp(arg[iarg],"linear") == 0) {
53 if (iarg+4 > narg) error->all(FLERR,"Illegal fix lb/momentum command");
54 linear = 1;
55 xflag = atoi(arg[iarg+1]);
56 yflag = atoi(arg[iarg+2]);
57 zflag = atoi(arg[iarg+3]);
58 iarg += 4;
59 } else error->all(FLERR,"Illegal fix lb/momentum command");
60 }
61
62 if (linear == 0)
63 error->all(FLERR,"Illegal fix lb/momentum command");
64
65 if (linear)
66 if (xflag < 0 || xflag > 1 || yflag < 0 || yflag > 1 ||
67 zflag < 0 || zflag > 1) error->all(FLERR,"Illegal fix lb/momentum command");
68
69 // cannot have 0 atoms in group
70
71 if (group->count(igroup) == 0.0)
72 error->all(FLERR,"Fix lb/momentum group has no atoms");
73
74 for (int ifix=0; ifix<modify->nfix; ifix++)
75 if (strcmp(modify->fix[ifix]->style,"lb/fluid")==0)
76 fix_lb_fluid = (FixLbFluid *)modify->fix[ifix];
77
78
79
80 }
81
82 /* ---------------------------------------------------------------------- */
83
setmask()84 int FixLbMomentum::setmask()
85 {
86 int mask = 0;
87 mask |= END_OF_STEP;
88 return mask;
89 }
90
91 /* ---------------------------------------------------------------------- */
92
init()93 void FixLbMomentum::init()
94 {
95 masstotal = group->mass(igroup);
96 }
97
98 /* ---------------------------------------------------------------------- */
99
end_of_step()100 void FixLbMomentum::end_of_step()
101 {
102 double masslb,masslbloc;
103 double momentumlbloc[3],momentumlb[3];
104 double vcmtotal[3];
105 const int numvel = fix_lb_fluid->numvel;
106 double etacov[19]; // = double etacov[numvel]; i.e. 15 or 19
107 double rho;
108
109 if (linear) {
110 double vcm[3];
111 group->vcm(igroup,masstotal,vcm);
112
113 // adjust velocities by vcm to zero linear momentum
114 // only adjust a component if flag is set
115
116 double **v = atom->v;
117 int *mask = atom->mask;
118 int nlocal = atom->nlocal;
119
120 int typeLB = fix_lb_fluid->typeLB;
121 int subNbx = fix_lb_fluid->subNbx;
122 int subNby = fix_lb_fluid->subNby;
123 int subNbz = fix_lb_fluid->subNbz;
124 double ***density_lb = fix_lb_fluid->density_lb;
125 double ****f_lb = fix_lb_fluid->f_lb;
126 double ****u_lb = fix_lb_fluid->u_lb;
127 double dm_lb = fix_lb_fluid->dm_lb;
128 double dx_lb = fix_lb_fluid->dx_lb;
129 double dt_lb = fix_lb_fluid->dt_lb;
130 double *Ng_lb = fix_lb_fluid->Ng_lb;
131 double *w_lb = fix_lb_fluid->w_lb;
132 double **mg_lb = fix_lb_fluid->mg_lb;
133 double ucmx,ucmy,ucmz;
134
135 //Calculate the total fluid mass and momentum.
136 masslbloc = 0.0;
137 momentumlbloc[0] = momentumlbloc[1] = momentumlbloc[2] = 0.0;
138
139 for (int i = 1; i<subNbx-1; i++)
140 for (int j = 1; j<subNby-1; j++)
141 for (int k = 1; k<subNbz-1; k++) {
142 masslbloc += density_lb[i][j][k];
143 momentumlbloc[0] += density_lb[i][j][k]*u_lb[i][j][k][0];
144 momentumlbloc[1] += density_lb[i][j][k]*u_lb[i][j][k][1];
145 momentumlbloc[2] += density_lb[i][j][k]*u_lb[i][j][k][2];
146 }
147 MPI_Allreduce(&masslbloc,&masslb,1,MPI_DOUBLE,MPI_SUM,world);
148 MPI_Allreduce(&momentumlbloc[0],&momentumlb[0],3,MPI_DOUBLE,MPI_SUM,world);
149
150 momentumlb[0] *= dm_lb*dx_lb/dt_lb;
151 momentumlb[1] *= dm_lb*dx_lb/dt_lb;
152 momentumlb[2] *= dm_lb*dx_lb/dt_lb;
153 masslb *= dm_lb;
154
155 //Calculate the center of mass velocity of the entire system.
156 vcmtotal[0] = (masstotal*vcm[0] + momentumlb[0])/(masslb + masstotal);
157 vcmtotal[1] = (masstotal*vcm[1] + momentumlb[1])/(masslb + masstotal);
158 vcmtotal[2] = (masstotal*vcm[2] + momentumlb[2])/(masslb + masstotal);
159
160 //Subtract vcm from the particles.
161 for (int i = 0; i < nlocal; i++)
162 if (mask[i] & groupbit) {
163 if (xflag) v[i][0] -= vcmtotal[0];
164 if (yflag) v[i][1] -= vcmtotal[1];
165 if (zflag) v[i][2] -= vcmtotal[2];
166 }
167
168 vcmtotal[0] *= dt_lb/dx_lb;
169 vcmtotal[1] *= dt_lb/dx_lb;
170 vcmtotal[2] *= dt_lb/dx_lb;
171
172 ucmx = ucmy = ucmz = 0.0;
173 double density_old;
174 double u_old[3];
175 int **e = fix_lb_fluid->e;
176
177 //Subtract vcm from the fluid.
178 for (int i=0; i<subNbx; i++)
179 for (int j=0; j<subNby; j++)
180 for (int k=0; k<subNbz; k++) {
181 rho = density_lb[i][j][k];
182 if (xflag) ucmx = vcmtotal[0];
183 if (yflag) ucmy = vcmtotal[1];
184 if (zflag) ucmz = vcmtotal[2];
185 if (numvel==15) {
186 etacov[0]=0.0;
187 etacov[1]=rho*ucmx;
188 etacov[2]=rho*ucmy;
189 etacov[3]=rho*ucmz;
190 etacov[4]=rho*(2.*u_lb[i][j][k][0]*ucmx-ucmx*ucmx);
191 etacov[5]=rho*(2.*u_lb[i][j][k][1]*ucmy-ucmy*ucmy);
192 etacov[6]=rho*(2.*u_lb[i][j][k][2]*ucmz-ucmz*ucmz);
193 etacov[7]=rho*(u_lb[i][j][k][0]*ucmy+u_lb[i][j][k][1]*ucmx-ucmx*ucmy);
194 etacov[8]=rho*(u_lb[i][j][k][1]*ucmz+u_lb[i][j][k][2]*ucmy-ucmy*ucmz);
195 etacov[9]=rho*(u_lb[i][j][k][0]*ucmz+u_lb[i][j][k][2]*ucmx-ucmx*ucmz);
196 etacov[10]=0.0;
197 etacov[11]=0.0;
198 etacov[12]=0.0;
199 etacov[13]=rho*(u_lb[i][j][k][0]*u_lb[i][j][k][1]*ucmz+u_lb[i][j][k][0]*ucmy*u_lb[i][j][k][2]-
200 u_lb[i][j][k][0]*ucmy*ucmz+ucmx*u_lb[i][j][k][1]*u_lb[i][j][k][2]-
201 ucmx*u_lb[i][j][k][1]*ucmz-ucmx*ucmy*u_lb[i][j][k][2]+
202 ucmx*ucmy*ucmz);
203 etacov[14]=0.0;
204 } else {
205 etacov[0] = 0.0;
206 etacov[1] = rho*ucmx;
207 etacov[2] = rho*ucmy;
208 etacov[3] = rho*ucmz;
209 etacov[4]=rho*(2.*u_lb[i][j][k][0]*ucmx-ucmx*ucmx);
210 etacov[5]=rho*(2.*u_lb[i][j][k][1]*ucmy-ucmy*ucmy);
211 etacov[6]=rho*(2.*u_lb[i][j][k][2]*ucmz-ucmz*ucmz);
212 etacov[7]=rho*(u_lb[i][j][k][0]*ucmy+u_lb[i][j][k][1]*ucmx-ucmx*ucmy);
213 etacov[8]=rho*(u_lb[i][j][k][0]*ucmz+u_lb[i][j][k][2]*ucmx-ucmx*ucmz);
214 etacov[9]=rho*(u_lb[i][j][k][1]*ucmz+u_lb[i][j][k][2]*ucmy-ucmy*ucmz);
215 etacov[10] = 0.0;
216 etacov[11] = 0.0;
217 etacov[12] = 0.0;
218 etacov[13] = 0.0;
219 etacov[14] = 0.0;
220 etacov[15] = 0.0;
221 etacov[16] = 0.0;
222 etacov[17] = 0.0;
223 etacov[18] = 0.0;
224 }
225
226 for (int l=0; l<numvel; l++)
227 for (int ii=0; ii<numvel; ii++) {
228 f_lb[i][j][k][l] -= w_lb[l]*mg_lb[ii][l]*etacov[ii]*Ng_lb[ii];
229 }
230
231
232 if (typeLB == 2) {
233 double ****feqold = fix_lb_fluid->feqold;
234 double ****feqoldn = fix_lb_fluid->feqoldn;
235 density_old = 0.0;
236 u_old[0] = u_old[1] = u_old[2] = 0.0;
237 for (int l=0; l<numvel; l++) {
238 density_old += feqold[i][j][k][l];
239 u_old[0] += feqold[i][j][k][l]*e[l][0];
240 u_old[1] += feqold[i][j][k][l]*e[l][1];
241 u_old[2] += feqold[i][j][k][l]*e[l][2];
242 }
243 u_old[0] = u_old[0]/density_old;
244 u_old[1] = u_old[1]/density_old;
245 u_old[2] = u_old[2]/density_old;
246
247 if (numvel==15) {
248 etacov[0]=0.0;
249 etacov[1]=density_old*ucmx;
250 etacov[2]=density_old*ucmy;
251 etacov[3]=density_old*ucmz;
252 etacov[4]=density_old*(2.*u_old[0]*ucmx-ucmx*ucmx);
253 etacov[5]=density_old*(2.*u_old[1]*ucmy-ucmy*ucmy);
254 etacov[6]=density_old*(2.*u_old[2]*ucmz-ucmz*ucmz);
255 etacov[7]=density_old*(u_old[0]*ucmy+u_old[1]*ucmx-ucmx*ucmy);
256 etacov[8]=density_old*(u_old[1]*ucmz+u_old[2]*ucmy-ucmy*ucmz);
257 etacov[9]=density_old*(u_old[0]*ucmz+u_old[2]*ucmx-ucmx*ucmz);
258 etacov[10]=0.0;
259 etacov[11]=0.0;
260 etacov[12]=0.0;
261 etacov[13]=density_old*(u_old[0]*u_old[1]*ucmz+u_old[0]*ucmy*u_old[2]-
262 u_old[0]*ucmy*ucmz+ucmx*u_old[1]*u_old[2]-
263 ucmx*u_old[1]*ucmz-ucmx*ucmy*u_old[2]+
264 ucmx*ucmy*ucmz);
265 etacov[14]=0.0;
266 } else {
267 etacov[0] = 0.0;
268 etacov[1] = density_old*ucmx;
269 etacov[2] = density_old*ucmy;
270 etacov[3] = density_old*ucmz;
271 etacov[4] = density_old*(2.*u_lb[i][j][k][0]*ucmx-ucmx*ucmx);
272 etacov[5] = density_old*(2.*u_lb[i][j][k][1]*ucmy-ucmy*ucmy);
273 etacov[6] = density_old*(2.*u_lb[i][j][k][2]*ucmz-ucmz*ucmz);
274 etacov[7] = density_old*(u_lb[i][j][k][0]*ucmy+u_lb[i][j][k][1]*ucmx-ucmx*ucmy);
275 etacov[8] = density_old*(u_lb[i][j][k][0]*ucmz+u_lb[i][j][k][2]*ucmx-ucmx*ucmz);
276 etacov[9] = density_old*(u_lb[i][j][k][1]*ucmz+u_lb[i][j][k][2]*ucmy-ucmy*ucmz);
277 etacov[10] = 0.0;
278 etacov[11] = 0.0;
279 etacov[12] = 0.0;
280 etacov[13] = 0.0;
281 etacov[14] = 0.0;
282 etacov[15] = 0.0;
283 etacov[16] = 0.0;
284 etacov[17] = 0.0;
285 etacov[18] = 0.0;
286 }
287
288 for (int l=0; l<numvel; l++)
289 for (int ii=0; ii<numvel; ii++) {
290 feqold[i][j][k][l] -= w_lb[l]*mg_lb[ii][l]*etacov[ii]*Ng_lb[ii];
291 feqoldn[i][j][k][l] -= w_lb[l]*mg_lb[ii][l]*etacov[ii]*Ng_lb[ii];
292 }
293 }
294 }
295 }
296
297 fix_lb_fluid->parametercalc_full();
298
299 }
300