1
/* ----------------------------------------------------------------------
2
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
3
http://lammps.sandia.gov, Sandia National Laboratories
4
Steve Plimpton, sjplimp@sandia.gov
6
Copyright (2003) Sandia Corporation. Under the terms of Contract
7
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
8
certain rights in this software. This software is distributed under
9
the GNU General Public License.
11
See the README file in the top-level LAMMPS directory.
12
------------------------------------------------------------------------- */
15
#include "verlet_kokkos.h"
20
#include "atom_kokkos.h"
21
#include "atom_masks.h"
40
using namespace LAMMPS_NS;
42
/* ---------------------------------------------------------------------- */
44
VerletKokkos::VerletKokkos(LAMMPS *lmp, int narg, char **arg) :
45
Verlet(lmp, narg, arg)
47
atomKK = (AtomKokkos *) atom;
50
/* ----------------------------------------------------------------------
52
------------------------------------------------------------------------- */
54
void VerletKokkos::setup()
57
if (comm->me == 0 && screen) fprintf(screen,"Setting up run ...\n");
58
update->setupflag = 1;
60
// setup domain, communication and neighboring
62
// build neighbor lists
64
atomKK->modified(Host,ALL_MASK);
67
modify->setup_pre_exchange();
69
atomKK->sync(Host,ALL_MASK);
70
atomKK->modified(Host,ALL_MASK);
71
if (triclinic) domain->x2lamda(atomKK->nlocal);
74
atomKK->sync(Host,ALL_MASK);
79
if (neighbor->style) neighbor->setup_bins();
83
if (atomKK->sortfreq > 0) atomKK->sort();
87
if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
89
atomKK->sync(Host,ALL_MASK);
91
domain->image_check();
92
domain->box_too_small_check();
93
modify->setup_pre_neighbor();
95
atomKK->modified(Host,ALL_MASK);
100
// compute all forces
102
ev_set(update->ntimestep);
104
modify->setup_pre_force(vflag);
106
if (pair_compute_flag) {
107
atomKK->sync(force->pair->execution_space,force->pair->datamask_read);
108
atomKK->modified(force->pair->execution_space,force->pair->datamask_modify);
109
force->pair->compute(eflag,vflag);
110
timer->stamp(TIME_PAIR);
112
else if (force->pair) force->pair->compute_dummy(eflag,vflag);
115
if (atomKK->molecular) {
117
atomKK->sync(force->bond->execution_space,force->bond->datamask_read);
118
atomKK->modified(force->bond->execution_space,force->bond->datamask_modify);
119
force->bond->compute(eflag,vflag);
122
atomKK->sync(force->angle->execution_space,force->angle->datamask_read);
123
atomKK->modified(force->angle->execution_space,force->angle->datamask_modify);
124
force->angle->compute(eflag,vflag);
126
if (force->dihedral) {
127
atomKK->sync(force->dihedral->execution_space,force->dihedral->datamask_read);
128
atomKK->modified(force->dihedral->execution_space,force->dihedral->datamask_modify);
129
force->dihedral->compute(eflag,vflag);
131
if (force->improper) {
132
atomKK->sync(force->improper->execution_space,force->improper->datamask_read);
133
atomKK->modified(force->improper->execution_space,force->improper->datamask_modify);
134
force->improper->compute(eflag,vflag);
136
timer->stamp(TIME_BOND);
140
force->kspace->setup();
141
if (kspace_compute_flag) {
142
atomKK->sync(force->kspace->execution_space,force->kspace->datamask_read);
143
atomKK->modified(force->kspace->execution_space,force->kspace->datamask_modify);
144
force->kspace->compute(eflag,vflag);
145
timer->stamp(TIME_KSPACE);
146
} else force->kspace->compute_dummy(eflag,vflag);
149
if (force->newton) comm->reverse_comm();
151
modify->setup(vflag);
153
update->setupflag = 0;
156
/* ----------------------------------------------------------------------
158
flag = 0 = just force calculation
159
flag = 1 = reneighbor and force calculation
160
------------------------------------------------------------------------- */
162
void VerletKokkos::setup_minimal(int flag)
164
update->setupflag = 1;
166
// setup domain, communication and neighboring
168
// build neighbor lists
171
atomKK->modified(Host,ALL_MASK);
173
modify->setup_pre_exchange();
175
atomKK->sync(Host,ALL_MASK);
176
atomKK->modified(Host,ALL_MASK);
178
if (triclinic) domain->x2lamda(atomKK->nlocal);
181
atomKK->sync(Host,ALL_MASK);
185
if (neighbor->style) neighbor->setup_bins();
188
if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
190
atomKK->sync(Host,ALL_MASK);
192
domain->image_check();
193
domain->box_too_small_check();
194
modify->setup_pre_neighbor();
196
atomKK->modified(Host,ALL_MASK);
199
neighbor->ncalls = 0;
202
// compute all forces
204
ev_set(update->ntimestep);
206
modify->setup_pre_force(vflag);
208
if (pair_compute_flag) {
209
atomKK->sync(force->pair->execution_space,force->pair->datamask_read);
210
atomKK->modified(force->pair->execution_space,force->pair->datamask_modify);
211
force->pair->compute(eflag,vflag);
212
timer->stamp(TIME_PAIR);
214
else if (force->pair) force->pair->compute_dummy(eflag,vflag);
217
if (atomKK->molecular) {
219
atomKK->sync(force->bond->execution_space,force->bond->datamask_read);
220
atomKK->modified(force->bond->execution_space,force->bond->datamask_modify);
221
force->bond->compute(eflag,vflag);
224
atomKK->sync(force->angle->execution_space,force->angle->datamask_read);
225
atomKK->modified(force->angle->execution_space,force->angle->datamask_modify);
226
force->angle->compute(eflag,vflag);
228
if (force->dihedral) {
229
atomKK->sync(force->dihedral->execution_space,force->dihedral->datamask_read);
230
atomKK->modified(force->dihedral->execution_space,force->dihedral->datamask_modify);
231
force->dihedral->compute(eflag,vflag);
233
if (force->improper) {
234
atomKK->sync(force->improper->execution_space,force->improper->datamask_read);
235
atomKK->modified(force->improper->execution_space,force->improper->datamask_modify);
236
force->improper->compute(eflag,vflag);
238
timer->stamp(TIME_BOND);
242
force->kspace->setup();
243
if (kspace_compute_flag) {
244
atomKK->sync(force->kspace->execution_space,force->kspace->datamask_read);
245
atomKK->modified(force->kspace->execution_space,force->kspace->datamask_modify);
246
force->kspace->compute(eflag,vflag);
247
timer->stamp(TIME_KSPACE);
248
} else force->kspace->compute_dummy(eflag,vflag);
251
if (force->newton) comm->reverse_comm();
253
modify->setup(vflag);
254
update->setupflag = 0;
257
/* ----------------------------------------------------------------------
259
------------------------------------------------------------------------- */
261
void VerletKokkos::run(int n)
266
int n_post_integrate = modify->n_post_integrate;
267
int n_pre_exchange = modify->n_pre_exchange;
268
int n_pre_neighbor = modify->n_pre_neighbor;
269
int n_pre_force = modify->n_pre_force;
270
int n_post_force = modify->n_post_force;
271
int n_end_of_step = modify->n_end_of_step;
273
if (atomKK->sortfreq > 0) sortflag = 1;
276
static double time = 0.0;
277
static int count = 0;
278
atomKK->sync(Device,ALL_MASK);
279
Kokkos::Impl::Timer ktimer;
281
for (int i = 0; i < n; i++) {
283
ntimestep = ++update->ntimestep;
286
// initial time integration
289
modify->initial_integrate(vflag);
290
time += ktimer.seconds();
291
if (n_post_integrate) modify->post_integrate();
293
// regular communication vs neighbor list rebuild
295
nflag = neighbor->decide();
299
comm->forward_comm();
300
timer->stamp(TIME_COMM);
303
//atomKK->sync(Host,ALL_MASK);
304
//atomKK->modified(Host,ALL_MASK);
306
if (n_pre_exchange) modify->pre_exchange();
308
//atomKK->sync(Host,ALL_MASK);
309
//atomKK->modified(Host,ALL_MASK);
310
if (triclinic) domain->x2lamda(atomKK->nlocal);
312
if (domain->box_change) {
315
if (neighbor->style) neighbor->setup_bins();
320
//atomKK->sync(Device,ALL_MASK);
321
//atomKK->modified(Device,ALL_MASK);
324
if (sortflag && ntimestep >= atomKK->nextsort) atomKK->sort();
328
//atomKK->sync(Host,ALL_MASK);
329
//atomKK->modified(Host,ALL_MASK);
331
if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
333
timer->stamp(TIME_COMM);
334
if (n_pre_neighbor) modify->pre_neighbor();
336
timer->stamp(TIME_NEIGHBOR);
339
// force computations
340
// important for pair to come before bonded contributions
341
// since some bonded potentials tally pairwise energy/virial
342
// and Pair:ev_tally() needs to be called before any tallying
346
//atomKK->k_x.sync<LMPHostType>();
347
//atomKK->k_f.sync<LMPHostType>();
348
//atomKK->k_f.modify<LMPHostType>();
349
if (n_pre_force) modify->pre_force(vflag);
353
if (pair_compute_flag) {
354
atomKK->sync(force->pair->execution_space,force->pair->datamask_read);
355
atomKK->modified(force->pair->execution_space,force->pair->datamask_modify);
356
force->pair->compute(eflag,vflag);
357
timer->stamp(TIME_PAIR);
360
if (atomKK->molecular) {
362
atomKK->sync(force->bond->execution_space,force->bond->datamask_read);
363
atomKK->modified(force->bond->execution_space,force->bond->datamask_modify);
364
force->bond->compute(eflag,vflag);
367
atomKK->sync(force->angle->execution_space,force->angle->datamask_read);
368
atomKK->modified(force->angle->execution_space,force->angle->datamask_modify);
369
force->angle->compute(eflag,vflag);
371
if (force->dihedral) {
372
atomKK->sync(force->dihedral->execution_space,force->dihedral->datamask_read);
373
atomKK->modified(force->dihedral->execution_space,force->dihedral->datamask_modify);
374
force->dihedral->compute(eflag,vflag);
376
if (force->improper) {
377
atomKK->sync(force->improper->execution_space,force->improper->datamask_read);
378
atomKK->modified(force->improper->execution_space,force->improper->datamask_modify);
379
force->improper->compute(eflag,vflag);
381
timer->stamp(TIME_BOND);
384
if (kspace_compute_flag) {
385
atomKK->sync(force->kspace->execution_space,force->kspace->datamask_read);
386
atomKK->modified(force->kspace->execution_space,force->kspace->datamask_modify);
387
force->kspace->compute(eflag,vflag);
388
timer->stamp(TIME_KSPACE);
391
// reverse communication of forces
393
if (force->newton) comm->reverse_comm();
394
timer->stamp(TIME_COMM);
396
// force modifications, final time integration, diagnostics
400
if (n_post_force) modify->post_force(vflag);
401
modify->final_integrate();
402
if (n_end_of_step) modify->end_of_step();
404
time += ktimer.seconds();
408
if (ntimestep == output->next) {
409
atomKK->sync(Host,ALL_MASK);
412
output->write(ntimestep);
413
timer->stamp(TIME_OUTPUT);
418
/* ----------------------------------------------------------------------
419
clear force on own & ghost atoms
420
clear other arrays as needed
421
------------------------------------------------------------------------- */
423
void VerletKokkos::force_clear()
427
if (external_force_clear) return;
429
// clear force on all particles
430
// if either newton flag is set, also include ghosts
431
// when using threads always clear all forces.
433
if (neighbor->includegroup == 0) {
435
if (force->newton) nall = atomKK->nlocal + atomKK->nghost;
436
else nall = atomKK->nlocal;
438
size_t nbytes = sizeof(double) * nall;
441
if (atomKK->k_f.modified_host > atomKK->k_f.modified_device) {
442
memset_kokkos(atomKK->k_f.view<LMPHostType>());
443
atomKK->modified(Host,F_MASK);
445
memset_kokkos(atomKK->k_f.view<LMPDeviceType>());
446
atomKK->modified(Device,F_MASK);
448
if (torqueflag) memset(&(atomKK->torque[0][0]),0,3*nbytes);
452
// neighbor includegroup flag is set
453
// clear force only on initial nfirst particles
454
// if either newton flag is set, also include ghosts
457
int nall = atomKK->nfirst;
458
if (atomKK->k_f.modified_host > atomKK->k_f.modified_device) {
459
memset_kokkos(atomKK->k_f.view<LMPHostType>());
460
atomKK->modified(Host,F_MASK);
462
memset_kokkos(atomKK->k_f.view<LMPDeviceType>());
463
atomKK->modified(Device,F_MASK);
466
double **torque = atomKK->torque;
467
for (i = 0; i < nall; i++) {
475
nall = atomKK->nlocal + atomKK->nghost;
478
double **torque = atomKK->torque;
479
for (i = atomKK->nlocal; i < nall; i++) {