~ubuntu-branches/debian/sid/lammps/sid

« back to all changes in this revision

Viewing changes to src/KOKKOS/verlet_kokkos.cpp

  • Committer: Package Import Robot
  • Author(s): Anton Gladky
  • Date: 2015-04-29 23:44:49 UTC
  • mfrom: (5.1.3 experimental)
  • Revision ID: package-import@ubuntu.com-20150429234449-mbhy9utku6hp6oq8
Tags: 0~20150313.gitfa668e1-1
Upload into unstable.

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
/* ----------------------------------------------------------------------
 
2
   LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
 
3
   http://lammps.sandia.gov, Sandia National Laboratories
 
4
   Steve Plimpton, sjplimp@sandia.gov
 
5
 
 
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.
 
10
 
 
11
   See the README file in the top-level LAMMPS directory.
 
12
------------------------------------------------------------------------- */
 
13
 
 
14
#include "string.h"
 
15
#include "verlet_kokkos.h"
 
16
#include "neighbor.h"
 
17
#include "domain.h"
 
18
#include "comm.h"
 
19
#include "atom.h"
 
20
#include "atom_kokkos.h"
 
21
#include "atom_masks.h"
 
22
#include "force.h"
 
23
#include "pair.h"
 
24
#include "bond.h"
 
25
#include "angle.h"
 
26
#include "dihedral.h"
 
27
#include "improper.h"
 
28
#include "kspace.h"
 
29
#include "output.h"
 
30
#include "update.h"
 
31
#include "modify.h"
 
32
#include "compute.h"
 
33
#include "fix.h"
 
34
#include "timer.h"
 
35
#include "memory.h"
 
36
#include "error.h"
 
37
 
 
38
#include <ctime>
 
39
 
 
40
using namespace LAMMPS_NS;
 
41
 
 
42
/* ---------------------------------------------------------------------- */
 
43
 
 
44
VerletKokkos::VerletKokkos(LAMMPS *lmp, int narg, char **arg) :
 
45
  Verlet(lmp, narg, arg) 
 
46
{
 
47
  atomKK = (AtomKokkos *) atom;
 
48
}
 
49
 
 
50
/* ----------------------------------------------------------------------
 
51
   setup before run
 
52
------------------------------------------------------------------------- */
 
53
 
 
54
void VerletKokkos::setup()
 
55
{
 
56
 
 
57
  if (comm->me == 0 && screen) fprintf(screen,"Setting up run ...\n");
 
58
  update->setupflag = 1;
 
59
 
 
60
  // setup domain, communication and neighboring
 
61
  // acquire ghosts
 
62
  // build neighbor lists
 
63
 
 
64
  atomKK->modified(Host,ALL_MASK);
 
65
 
 
66
  atomKK->setup();
 
67
  modify->setup_pre_exchange();
 
68
      // debug
 
69
  atomKK->sync(Host,ALL_MASK);
 
70
  atomKK->modified(Host,ALL_MASK);
 
71
  if (triclinic) domain->x2lamda(atomKK->nlocal);
 
72
  domain->pbc();
 
73
 
 
74
  atomKK->sync(Host,ALL_MASK);
 
75
 
 
76
 
 
77
  domain->reset_box();
 
78
  comm->setup();
 
79
  if (neighbor->style) neighbor->setup_bins();
 
80
 
 
81
  comm->exchange();
 
82
 
 
83
  if (atomKK->sortfreq > 0) atomKK->sort();
 
84
 
 
85
  comm->borders();
 
86
 
 
87
  if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
 
88
 
 
89
  atomKK->sync(Host,ALL_MASK);
 
90
 
 
91
  domain->image_check();
 
92
  domain->box_too_small_check();
 
93
  modify->setup_pre_neighbor();
 
94
 
 
95
  atomKK->modified(Host,ALL_MASK);
 
96
 
 
97
  neighbor->build();
 
98
  neighbor->ncalls = 0;
 
99
 
 
100
  // compute all forces
 
101
 
 
102
  ev_set(update->ntimestep);
 
103
  force_clear();
 
104
  modify->setup_pre_force(vflag);
 
105
 
 
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);
 
111
  }
 
112
  else if (force->pair) force->pair->compute_dummy(eflag,vflag);
 
113
 
 
114
 
 
115
  if (atomKK->molecular) {
 
116
    if (force->bond) {
 
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);
 
120
    }
 
121
    if (force->angle) {
 
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);
 
125
    }
 
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);
 
130
    }
 
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);
 
135
    }
 
136
    timer->stamp(TIME_BOND);
 
137
  }
 
138
 
 
139
  if(force->kspace) {
 
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);
 
147
  }
 
148
 
 
149
  if (force->newton) comm->reverse_comm();
 
150
 
 
151
  modify->setup(vflag);
 
152
  output->setup();
 
153
  update->setupflag = 0;
 
154
}
 
155
 
 
156
/* ----------------------------------------------------------------------
 
157
   setup without output
 
158
   flag = 0 = just force calculation
 
159
   flag = 1 = reneighbor and force calculation
 
160
------------------------------------------------------------------------- */
 
161
 
 
162
void VerletKokkos::setup_minimal(int flag)
 
163
{
 
164
  update->setupflag = 1;
 
165
 
 
166
  // setup domain, communication and neighboring
 
167
  // acquire ghosts
 
168
  // build neighbor lists
 
169
 
 
170
  if (flag) {
 
171
    atomKK->modified(Host,ALL_MASK);
 
172
 
 
173
    modify->setup_pre_exchange();
 
174
      // debug
 
175
      atomKK->sync(Host,ALL_MASK);
 
176
      atomKK->modified(Host,ALL_MASK);
 
177
 
 
178
    if (triclinic) domain->x2lamda(atomKK->nlocal);
 
179
    domain->pbc();
 
180
 
 
181
    atomKK->sync(Host,ALL_MASK);
 
182
 
 
183
    domain->reset_box();
 
184
    comm->setup();
 
185
    if (neighbor->style) neighbor->setup_bins();
 
186
    comm->exchange();
 
187
    comm->borders();
 
188
    if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
 
189
 
 
190
    atomKK->sync(Host,ALL_MASK);
 
191
 
 
192
    domain->image_check();
 
193
    domain->box_too_small_check();
 
194
    modify->setup_pre_neighbor();
 
195
 
 
196
    atomKK->modified(Host,ALL_MASK);
 
197
 
 
198
    neighbor->build();
 
199
    neighbor->ncalls = 0;
 
200
  }
 
201
 
 
202
  // compute all forces
 
203
 
 
204
  ev_set(update->ntimestep);
 
205
  force_clear();
 
206
  modify->setup_pre_force(vflag);
 
207
 
 
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);
 
213
  }
 
214
  else if (force->pair) force->pair->compute_dummy(eflag,vflag);
 
215
 
 
216
 
 
217
  if (atomKK->molecular) {
 
218
    if (force->bond) {
 
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);
 
222
    }
 
223
    if (force->angle) {
 
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);
 
227
    }
 
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);
 
232
    }
 
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);
 
237
    }
 
238
    timer->stamp(TIME_BOND);
 
239
  }
 
240
 
 
241
  if(force->kspace) {
 
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);
 
249
  }
 
250
 
 
251
  if (force->newton) comm->reverse_comm();
 
252
 
 
253
  modify->setup(vflag);
 
254
  update->setupflag = 0;
 
255
}
 
256
 
 
257
/* ----------------------------------------------------------------------
 
258
   run for N steps
 
259
------------------------------------------------------------------------- */
 
260
 
 
261
void VerletKokkos::run(int n)
 
262
{
 
263
  bigint ntimestep;
 
264
  int nflag,sortflag;
 
265
 
 
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;
 
272
 
 
273
  if (atomKK->sortfreq > 0) sortflag = 1;
 
274
  else sortflag = 0;
 
275
 
 
276
  static double time = 0.0;
 
277
  static int count = 0;
 
278
  atomKK->sync(Device,ALL_MASK);
 
279
  Kokkos::Impl::Timer ktimer;
 
280
 
 
281
  for (int i = 0; i < n; i++) {
 
282
 
 
283
    ntimestep = ++update->ntimestep;
 
284
    ev_set(ntimestep);
 
285
 
 
286
    // initial time integration
 
287
 
 
288
    ktimer.reset();
 
289
    modify->initial_integrate(vflag);
 
290
    time += ktimer.seconds();
 
291
    if (n_post_integrate) modify->post_integrate();
 
292
 
 
293
    // regular communication vs neighbor list rebuild
 
294
 
 
295
    nflag = neighbor->decide();
 
296
 
 
297
    if (nflag == 0) {
 
298
      timer->stamp();
 
299
      comm->forward_comm();
 
300
      timer->stamp(TIME_COMM);
 
301
    } else {
 
302
      // added debug
 
303
      //atomKK->sync(Host,ALL_MASK);
 
304
      //atomKK->modified(Host,ALL_MASK);
 
305
 
 
306
      if (n_pre_exchange) modify->pre_exchange();
 
307
      // debug
 
308
      //atomKK->sync(Host,ALL_MASK);
 
309
      //atomKK->modified(Host,ALL_MASK);
 
310
      if (triclinic) domain->x2lamda(atomKK->nlocal);
 
311
      domain->pbc();
 
312
      if (domain->box_change) {
 
313
        domain->reset_box();
 
314
        comm->setup();
 
315
        if (neighbor->style) neighbor->setup_bins();
 
316
      }
 
317
      timer->stamp();
 
318
 
 
319
      // added debug
 
320
      //atomKK->sync(Device,ALL_MASK);
 
321
      //atomKK->modified(Device,ALL_MASK);
 
322
 
 
323
      comm->exchange();
 
324
      if (sortflag && ntimestep >= atomKK->nextsort) atomKK->sort();
 
325
      comm->borders();
 
326
 
 
327
      // added debug
 
328
      //atomKK->sync(Host,ALL_MASK);
 
329
      //atomKK->modified(Host,ALL_MASK);
 
330
 
 
331
      if (triclinic) domain->lamda2x(atomKK->nlocal+atomKK->nghost);
 
332
 
 
333
      timer->stamp(TIME_COMM);
 
334
      if (n_pre_neighbor) modify->pre_neighbor();
 
335
      neighbor->build();
 
336
      timer->stamp(TIME_NEIGHBOR);
 
337
    }
 
338
 
 
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
 
343
 
 
344
    force_clear();
 
345
    // added for debug
 
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);
 
350
 
 
351
    timer->stamp();
 
352
 
 
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);
 
358
    }
 
359
 
 
360
    if (atomKK->molecular) {
 
361
      if (force->bond) {
 
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);
 
365
      }
 
366
      if (force->angle) {
 
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);
 
370
      }
 
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);
 
375
      }
 
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);
 
380
      }
 
381
      timer->stamp(TIME_BOND);
 
382
    }
 
383
 
 
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);
 
389
    }
 
390
 
 
391
    // reverse communication of forces
 
392
 
 
393
    if (force->newton) comm->reverse_comm();
 
394
    timer->stamp(TIME_COMM);
 
395
 
 
396
    // force modifications, final time integration, diagnostics
 
397
 
 
398
    ktimer.reset();
 
399
 
 
400
    if (n_post_force) modify->post_force(vflag);
 
401
    modify->final_integrate();
 
402
    if (n_end_of_step) modify->end_of_step();
 
403
 
 
404
    time += ktimer.seconds();
 
405
 
 
406
    // all output
 
407
 
 
408
    if (ntimestep == output->next) {
 
409
       atomKK->sync(Host,ALL_MASK);
 
410
 
 
411
      timer->stamp();
 
412
      output->write(ntimestep);
 
413
      timer->stamp(TIME_OUTPUT);
 
414
    }
 
415
  }
 
416
}
 
417
 
 
418
/* ----------------------------------------------------------------------
 
419
   clear force on own & ghost atoms
 
420
   clear other arrays as needed
 
421
------------------------------------------------------------------------- */
 
422
 
 
423
void VerletKokkos::force_clear()
 
424
{
 
425
  int i;
 
426
 
 
427
  if (external_force_clear) return;
 
428
 
 
429
  // clear force on all particles
 
430
  // if either newton flag is set, also include ghosts
 
431
  // when using threads always clear all forces.
 
432
 
 
433
  if (neighbor->includegroup == 0) {
 
434
    int nall;
 
435
    if (force->newton) nall = atomKK->nlocal + atomKK->nghost;
 
436
    else nall = atomKK->nlocal;
 
437
 
 
438
    size_t nbytes = sizeof(double) * nall;
 
439
 
 
440
    if (nbytes) {
 
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);
 
444
      } else {
 
445
        memset_kokkos(atomKK->k_f.view<LMPDeviceType>());
 
446
        atomKK->modified(Device,F_MASK);
 
447
      }
 
448
      if (torqueflag)  memset(&(atomKK->torque[0][0]),0,3*nbytes);
 
449
 
 
450
    }
 
451
 
 
452
  // neighbor includegroup flag is set
 
453
  // clear force only on initial nfirst particles
 
454
  // if either newton flag is set, also include ghosts
 
455
 
 
456
  } else {
 
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);
 
461
    } else {
 
462
      memset_kokkos(atomKK->k_f.view<LMPDeviceType>());
 
463
      atomKK->modified(Device,F_MASK);
 
464
    }
 
465
    if (torqueflag) {
 
466
      double **torque = atomKK->torque;
 
467
      for (i = 0; i < nall; i++) {
 
468
        torque[i][0] = 0.0;
 
469
        torque[i][1] = 0.0;
 
470
        torque[i][2] = 0.0;
 
471
      }
 
472
    }
 
473
 
 
474
    if (force->newton) {
 
475
      nall = atomKK->nlocal + atomKK->nghost;
 
476
 
 
477
      if (torqueflag) {
 
478
        double **torque = atomKK->torque;
 
479
        for (i = atomKK->nlocal; i < nall; i++) {
 
480
          torque[i][0] = 0.0;
 
481
          torque[i][1] = 0.0;
 
482
          torque[i][2] = 0.0;
 
483
        }
 
484
      }
 
485
 
 
486
    }
 
487
  }
 
488
}