Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F88587810
pair_lj_charmm_coul_long_kokkos.cpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sat, Oct 19, 14:57
Size
16 KB
Mime Type
text/x-c
Expires
Mon, Oct 21, 14:57 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
21794071
Attached To
rLAMMPS lammps
pair_lj_charmm_coul_long_kokkos.cpp
View Options
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
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: Ray Shan (SNL)
------------------------------------------------------------------------- */
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "pair_lj_charmm_coul_long_kokkos.h"
#include "kokkos.h"
#include "atom_kokkos.h"
#include "comm.h"
#include "force.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "neigh_request.h"
#include "update.h"
#include "integrate.h"
#include "respa.h"
#include "math_const.h"
#include "memory.h"
#include "error.h"
#include "atom_masks.h"
using
namespace
LAMMPS_NS
;
using
namespace
MathConst
;
#define KOKKOS_CUDA_MAX_THREADS 256
#define KOKKOS_CUDA_MIN_BLOCKS 8
#define EWALD_F 1.12837917
#define EWALD_P 0.3275911
#define A1 0.254829592
#define A2 -0.284496736
#define A3 1.421413741
#define A4 -1.453152027
#define A5 1.061405429
/* ---------------------------------------------------------------------- */
template
<
class
DeviceType
>
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
PairLJCharmmCoulLongKokkos
(
LAMMPS
*
lmp
)
:
PairLJCharmmCoulLong
(
lmp
)
{
respa_enable
=
0
;
atomKK
=
(
AtomKokkos
*
)
atom
;
execution_space
=
ExecutionSpaceFromDevice
<
DeviceType
>::
space
;
datamask_read
=
X_MASK
|
F_MASK
|
TYPE_MASK
|
Q_MASK
|
ENERGY_MASK
|
VIRIAL_MASK
;
datamask_modify
=
F_MASK
|
ENERGY_MASK
|
VIRIAL_MASK
;
cutsq
=
NULL
;
cut_ljsq
=
0.0
;
cut_coulsq
=
0.0
;
}
/* ---------------------------------------------------------------------- */
template
<
class
DeviceType
>
PairLJCharmmCoulLongKokkos
<
DeviceType
>::~
PairLJCharmmCoulLongKokkos
()
{
if
(
!
copymode
)
{
memory
->
destroy_kokkos
(
k_eatom
,
eatom
);
memory
->
destroy_kokkos
(
k_vatom
,
vatom
);
k_cutsq
=
DAT
::
tdual_ffloat_2d
();
k_cut_ljsq
=
DAT
::
tdual_ffloat_2d
();
k_cut_coulsq
=
DAT
::
tdual_ffloat_2d
();
memory
->
sfree
(
cutsq
);
//memory->sfree(cut_ljsq);
//memory->sfree(cut_coulsq);
eatom
=
NULL
;
vatom
=
NULL
;
cutsq
=
NULL
;
cut_ljsq
=
0.0
;
cut_coulsq
=
0.0
;
}
}
/* ---------------------------------------------------------------------- */
template
<
class
DeviceType
>
void
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
cleanup_copy
()
{
// WHY needed: this prevents parent copy from deallocating any arrays
allocated
=
0
;
cutsq
=
NULL
;
cut_ljsq
=
0.0
;
eatom
=
NULL
;
vatom
=
NULL
;
ftable
=
NULL
;
}
/* ---------------------------------------------------------------------- */
template
<
class
DeviceType
>
void
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
compute
(
int
eflag_in
,
int
vflag_in
)
{
eflag
=
eflag_in
;
vflag
=
vflag_in
;
if
(
neighflag
==
FULL
||
neighflag
==
FULLCLUSTER
)
no_virial_fdotr_compute
=
1
;
if
(
eflag
||
vflag
)
ev_setup
(
eflag
,
vflag
);
else
evflag
=
vflag_fdotr
=
0
;
atomKK
->
sync
(
execution_space
,
datamask_read
);
k_cutsq
.
template
sync
<
DeviceType
>
();
k_cut_ljsq
.
template
sync
<
DeviceType
>
();
k_cut_coulsq
.
template
sync
<
DeviceType
>
();
k_params
.
template
sync
<
DeviceType
>
();
if
(
eflag
||
vflag
)
atomKK
->
modified
(
execution_space
,
datamask_modify
);
else
atomKK
->
modified
(
execution_space
,
F_MASK
);
x
=
atomKK
->
k_x
.
view
<
DeviceType
>
();
c_x
=
atomKK
->
k_x
.
view
<
DeviceType
>
();
f
=
atomKK
->
k_f
.
view
<
DeviceType
>
();
q
=
atomKK
->
k_q
.
view
<
DeviceType
>
();
type
=
atomKK
->
k_type
.
view
<
DeviceType
>
();
nlocal
=
atom
->
nlocal
;
nall
=
atom
->
nlocal
+
atom
->
nghost
;
special_lj
[
0
]
=
force
->
special_lj
[
0
];
special_lj
[
1
]
=
force
->
special_lj
[
1
];
special_lj
[
2
]
=
force
->
special_lj
[
2
];
special_lj
[
3
]
=
force
->
special_lj
[
3
];
special_coul
[
0
]
=
force
->
special_coul
[
0
];
special_coul
[
1
]
=
force
->
special_coul
[
1
];
special_coul
[
2
]
=
force
->
special_coul
[
2
];
special_coul
[
3
]
=
force
->
special_coul
[
3
];
qqrd2e
=
force
->
qqrd2e
;
newton_pair
=
force
->
newton_pair
;
// loop over neighbors of my atoms
copymode
=
1
;
EV_FLOAT
ev
;
if
(
ncoultablebits
)
ev
=
pair_compute
<
PairLJCharmmCoulLongKokkos
<
DeviceType
>
,
CoulLongTable
<
1
>
>
(
this
,(
NeighListKokkos
<
DeviceType
>*
)
list
);
else
ev
=
pair_compute
<
PairLJCharmmCoulLongKokkos
<
DeviceType
>
,
CoulLongTable
<
0
>
>
(
this
,(
NeighListKokkos
<
DeviceType
>*
)
list
);
DeviceType
::
fence
();
if
(
eflag
)
{
eng_vdwl
+=
ev
.
evdwl
;
eng_coul
+=
ev
.
ecoul
;
}
if
(
vflag_global
)
{
virial
[
0
]
+=
ev
.
v
[
0
];
virial
[
1
]
+=
ev
.
v
[
1
];
virial
[
2
]
+=
ev
.
v
[
2
];
virial
[
3
]
+=
ev
.
v
[
3
];
virial
[
4
]
+=
ev
.
v
[
4
];
virial
[
5
]
+=
ev
.
v
[
5
];
}
if
(
vflag_fdotr
)
pair_virial_fdotr_compute
(
this
);
copymode
=
0
;
}
/* ----------------------------------------------------------------------
compute LJ CHARMM pair force between atoms i and j
---------------------------------------------------------------------- */
template
<
class
DeviceType
>
template
<
bool
STACKPARAMS
,
class
Specialisation
>
KOKKOS_INLINE_FUNCTION
F_FLOAT
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
compute_fpair
(
const
F_FLOAT
&
rsq
,
const
int
&
i
,
const
int
&
j
,
const
int
&
itype
,
const
int
&
jtype
)
const
{
const
F_FLOAT
r2inv
=
1.0
/
rsq
;
const
F_FLOAT
r6inv
=
r2inv
*
r2inv
*
r2inv
;
F_FLOAT
forcelj
,
switch1
,
switch2
,
englj
;
forcelj
=
r6inv
*
((
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj1:
params
(
itype
,
jtype
).
lj1
)
*
r6inv
-
(
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj2:
params
(
itype
,
jtype
).
lj2
));
if
(
rsq
>
cut_lj_innersq
)
{
switch1
=
(
cut_ljsq
-
rsq
)
*
(
cut_ljsq
-
rsq
)
*
(
cut_ljsq
+
2.0
*
rsq
-
3.0
*
cut_lj_innersq
)
/
denom_lj
;
switch2
=
12.0
*
rsq
*
(
cut_ljsq
-
rsq
)
*
(
rsq
-
cut_lj_innersq
)
/
denom_lj
;
englj
=
r6inv
*
((
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj3:
params
(
itype
,
jtype
).
lj3
)
*
r6inv
-
(
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj4:
params
(
itype
,
jtype
).
lj4
));
forcelj
=
forcelj
*
switch1
+
englj
*
switch2
;
}
return
forcelj
*
r2inv
;
}
/* ----------------------------------------------------------------------
compute LJ CHARMM pair potential energy between atoms i and j
---------------------------------------------------------------------- */
template
<
class
DeviceType
>
template
<
bool
STACKPARAMS
,
class
Specialisation
>
KOKKOS_INLINE_FUNCTION
F_FLOAT
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
compute_evdwl
(
const
F_FLOAT
&
rsq
,
const
int
&
i
,
const
int
&
j
,
const
int
&
itype
,
const
int
&
jtype
)
const
{
const
F_FLOAT
r2inv
=
1.0
/
rsq
;
const
F_FLOAT
r6inv
=
r2inv
*
r2inv
*
r2inv
;
F_FLOAT
englj
,
switch1
;
englj
=
r6inv
*
((
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj3:
params
(
itype
,
jtype
).
lj3
)
*
r6inv
-
(
STACKPARAMS
?
m_params
[
itype
][
jtype
].
lj4:
params
(
itype
,
jtype
).
lj4
));
if
(
rsq
>
cut_lj_innersq
)
{
switch1
=
(
cut_ljsq
-
rsq
)
*
(
cut_ljsq
-
rsq
)
*
(
cut_ljsq
+
2.0
*
rsq
-
3.0
*
cut_lj_innersq
)
/
denom_lj
;
englj
*=
switch1
;
}
return
englj
;
}
/* ----------------------------------------------------------------------
compute coulomb pair force between atoms i and j
---------------------------------------------------------------------- */
template
<
class
DeviceType
>
template
<
bool
STACKPARAMS
,
class
Specialisation
>
KOKKOS_INLINE_FUNCTION
F_FLOAT
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
compute_fcoul
(
const
F_FLOAT
&
rsq
,
const
int
&
i
,
const
int
&
j
,
const
int
&
itype
,
const
int
&
jtype
,
const
F_FLOAT
&
factor_coul
,
const
F_FLOAT
&
qtmp
)
const
{
if
(
Specialisation
::
DoTable
&&
rsq
>
tabinnersq
)
{
union_int_float_t
rsq_lookup
;
rsq_lookup
.
f
=
rsq
;
const
int
itable
=
(
rsq_lookup
.
i
&
ncoulmask
)
>>
ncoulshiftbits
;
const
F_FLOAT
fraction
=
(
rsq_lookup
.
f
-
d_rtable
[
itable
])
*
d_drtable
[
itable
];
const
F_FLOAT
table
=
d_ftable
[
itable
]
+
fraction
*
d_dftable
[
itable
];
F_FLOAT
forcecoul
=
qtmp
*
q
[
j
]
*
table
;
if
(
factor_coul
<
1.0
)
{
const
F_FLOAT
table
=
d_ctable
[
itable
]
+
fraction
*
d_dctable
[
itable
];
const
F_FLOAT
prefactor
=
qtmp
*
q
[
j
]
*
table
;
forcecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
}
return
forcecoul
/
rsq
;
}
else
{
const
F_FLOAT
r
=
sqrt
(
rsq
);
const
F_FLOAT
grij
=
g_ewald
*
r
;
const
F_FLOAT
expm2
=
exp
(
-
grij
*
grij
);
const
F_FLOAT
t
=
1.0
/
(
1.0
+
EWALD_P
*
grij
);
const
F_FLOAT
rinv
=
1.0
/
r
;
const
F_FLOAT
erfc
=
t
*
(
A1
+
t
*
(
A2
+
t
*
(
A3
+
t
*
(
A4
+
t
*
A5
))))
*
expm2
;
const
F_FLOAT
prefactor
=
qqrd2e
*
qtmp
*
q
[
j
]
*
rinv
;
F_FLOAT
forcecoul
=
prefactor
*
(
erfc
+
EWALD_F
*
grij
*
expm2
);
if
(
factor_coul
<
1.0
)
forcecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
return
forcecoul
*
rinv
*
rinv
;
}
}
/* ----------------------------------------------------------------------
compute coulomb pair potential energy between atoms i and j
---------------------------------------------------------------------- */
template
<
class
DeviceType
>
template
<
bool
STACKPARAMS
,
class
Specialisation
>
KOKKOS_INLINE_FUNCTION
F_FLOAT
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
compute_ecoul
(
const
F_FLOAT
&
rsq
,
const
int
&
i
,
const
int
&
j
,
const
int
&
itype
,
const
int
&
jtype
,
const
F_FLOAT
&
factor_coul
,
const
F_FLOAT
&
qtmp
)
const
{
if
(
Specialisation
::
DoTable
&&
rsq
>
tabinnersq
)
{
union_int_float_t
rsq_lookup
;
rsq_lookup
.
f
=
rsq
;
const
int
itable
=
(
rsq_lookup
.
i
&
ncoulmask
)
>>
ncoulshiftbits
;
const
F_FLOAT
fraction
=
(
rsq_lookup
.
f
-
d_rtable
[
itable
])
*
d_drtable
[
itable
];
const
F_FLOAT
table
=
d_etable
[
itable
]
+
fraction
*
d_detable
[
itable
];
F_FLOAT
ecoul
=
qtmp
*
q
[
j
]
*
table
;
if
(
factor_coul
<
1.0
)
{
const
F_FLOAT
table
=
d_ctable
[
itable
]
+
fraction
*
d_dctable
[
itable
];
const
F_FLOAT
prefactor
=
qtmp
*
q
[
j
]
*
table
;
ecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
}
return
ecoul
;
}
else
{
const
F_FLOAT
r
=
sqrt
(
rsq
);
const
F_FLOAT
grij
=
g_ewald
*
r
;
const
F_FLOAT
expm2
=
exp
(
-
grij
*
grij
);
const
F_FLOAT
t
=
1.0
/
(
1.0
+
EWALD_P
*
grij
);
const
F_FLOAT
erfc
=
t
*
(
A1
+
t
*
(
A2
+
t
*
(
A3
+
t
*
(
A4
+
t
*
A5
))))
*
expm2
;
const
F_FLOAT
prefactor
=
qqrd2e
*
qtmp
*
q
[
j
]
/
r
;
F_FLOAT
ecoul
=
prefactor
*
erfc
;
if
(
factor_coul
<
1.0
)
ecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
return
ecoul
;
}
}
/* ----------------------------------------------------------------------
allocate all arrays
------------------------------------------------------------------------- */
template
<
class
DeviceType
>
void
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
allocate
()
{
PairLJCharmmCoulLong
::
allocate
();
int
n
=
atom
->
ntypes
;
memory
->
destroy
(
cutsq
);
memory
->
create_kokkos
(
k_cutsq
,
cutsq
,
n
+
1
,
n
+
1
,
"pair:cutsq"
);
d_cutsq
=
k_cutsq
.
template
view
<
DeviceType
>
();
//memory->destroy(cut_ljsq);
memory
->
create_kokkos
(
k_cut_ljsq
,
n
+
1
,
n
+
1
,
"pair:cut_ljsq"
);
d_cut_ljsq
=
k_cut_ljsq
.
template
view
<
DeviceType
>
();
memory
->
create_kokkos
(
k_cut_coulsq
,
n
+
1
,
n
+
1
,
"pair:cut_coulsq"
);
d_cut_coulsq
=
k_cut_coulsq
.
template
view
<
DeviceType
>
();
k_params
=
Kokkos
::
DualView
<
params_lj_coul
**
,
Kokkos
::
LayoutRight
,
DeviceType
>
(
"PairLJCharmmCoulLong::params"
,
n
+
1
,
n
+
1
);
params
=
k_params
.
d_view
;
}
template
<
class
DeviceType
>
void
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
init_tables
(
double
cut_coul
,
double
*
cut_respa
)
{
Pair
::
init_tables
(
cut_coul
,
cut_respa
);
typedef
typename
ArrayTypes
<
DeviceType
>::
t_ffloat_1d
table_type
;
typedef
typename
ArrayTypes
<
LMPHostType
>::
t_ffloat_1d
host_table_type
;
int
ntable
=
1
;
for
(
int
i
=
0
;
i
<
ncoultablebits
;
i
++
)
ntable
*=
2
;
// Copy rtable and drtable
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
rtable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_rtable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
drtable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_drtable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
// Copy ftable and dftable
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
ftable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_ftable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
dftable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_dftable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
// Copy ctable and dctable
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
ctable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_ctable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
dctable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_dctable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
// Copy etable and detable
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
etable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_etable
=
d_table
;
}
{
host_table_type
h_table
(
"HostTable"
,
ntable
);
table_type
d_table
(
"DeviceTable"
,
ntable
);
for
(
int
i
=
0
;
i
<
ntable
;
i
++
)
{
h_table
(
i
)
=
detable
[
i
];
}
Kokkos
::
deep_copy
(
d_table
,
h_table
);
d_detable
=
d_table
;
}
}
/* ----------------------------------------------------------------------
init specific to this pair style
------------------------------------------------------------------------- */
template
<
class
DeviceType
>
void
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
init_style
()
{
PairLJCharmmCoulLong
::
init_style
();
// error if rRESPA with inner levels
if
(
update
->
whichflag
==
1
&&
strstr
(
update
->
integrate_style
,
"respa"
))
{
int
respa
=
0
;
if
(((
Respa
*
)
update
->
integrate
)
->
level_inner
>=
0
)
respa
=
1
;
if
(((
Respa
*
)
update
->
integrate
)
->
level_middle
>=
0
)
respa
=
2
;
if
(
respa
)
error
->
all
(
FLERR
,
"Cannot use Kokkos pair style with rRESPA inner/middle"
);
}
// irequest = neigh request made by parent class
neighflag
=
lmp
->
kokkos
->
neighflag
;
int
irequest
=
neighbor
->
nrequest
-
1
;
neighbor
->
requests
[
irequest
]
->
kokkos_host
=
Kokkos
::
Impl
::
is_same
<
DeviceType
,
LMPHostType
>::
value
&&
!
Kokkos
::
Impl
::
is_same
<
DeviceType
,
LMPDeviceType
>::
value
;
neighbor
->
requests
[
irequest
]
->
kokkos_device
=
Kokkos
::
Impl
::
is_same
<
DeviceType
,
LMPDeviceType
>::
value
;
if
(
neighflag
==
FULL
)
{
neighbor
->
requests
[
irequest
]
->
full
=
1
;
neighbor
->
requests
[
irequest
]
->
half
=
0
;
neighbor
->
requests
[
irequest
]
->
full_cluster
=
0
;
}
else
if
(
neighflag
==
HALF
||
neighflag
==
HALFTHREAD
)
{
neighbor
->
requests
[
irequest
]
->
full
=
0
;
neighbor
->
requests
[
irequest
]
->
half
=
1
;
neighbor
->
requests
[
irequest
]
->
full_cluster
=
0
;
}
else
{
error
->
all
(
FLERR
,
"Cannot use chosen neighbor list style with lj/charmm/coul/long/kk"
);
}
}
/* ----------------------------------------------------------------------
init for one type pair i,j and corresponding j,i
------------------------------------------------------------------------- */
template
<
class
DeviceType
>
double
PairLJCharmmCoulLongKokkos
<
DeviceType
>::
init_one
(
int
i
,
int
j
)
{
double
cutone
=
PairLJCharmmCoulLong
::
init_one
(
i
,
j
);
double
cut_ljsqm
=
cut_ljsq
;
double
cut_coulsqm
=
cut_coulsq
;
k_params
.
h_view
(
i
,
j
).
lj1
=
lj1
[
i
][
j
];
k_params
.
h_view
(
i
,
j
).
lj2
=
lj2
[
i
][
j
];
k_params
.
h_view
(
i
,
j
).
lj3
=
lj3
[
i
][
j
];
k_params
.
h_view
(
i
,
j
).
lj4
=
lj4
[
i
][
j
];
//k_params.h_view(i,j).offset = offset[i][j];
k_params
.
h_view
(
i
,
j
).
cut_ljsq
=
cut_ljsqm
;
k_params
.
h_view
(
i
,
j
).
cut_coulsq
=
cut_coulsqm
;
k_params
.
h_view
(
j
,
i
)
=
k_params
.
h_view
(
i
,
j
);
if
(
i
<
MAX_TYPES_STACKPARAMS
+
1
&&
j
<
MAX_TYPES_STACKPARAMS
+
1
)
{
m_params
[
i
][
j
]
=
m_params
[
j
][
i
]
=
k_params
.
h_view
(
i
,
j
);
m_cutsq
[
j
][
i
]
=
m_cutsq
[
i
][
j
]
=
cutone
*
cutone
;
m_cut_ljsq
[
j
][
i
]
=
m_cut_ljsq
[
i
][
j
]
=
cut_ljsqm
;
m_cut_coulsq
[
j
][
i
]
=
m_cut_coulsq
[
i
][
j
]
=
cut_coulsqm
;
}
k_cutsq
.
h_view
(
i
,
j
)
=
cutone
*
cutone
;
k_cutsq
.
h_view
(
j
,
i
)
=
k_cutsq
.
h_view
(
i
,
j
);
k_cutsq
.
template
modify
<
LMPHostType
>
();
k_cut_ljsq
.
h_view
(
i
,
j
)
=
cut_ljsqm
;
k_cut_ljsq
.
h_view
(
j
,
i
)
=
k_cut_ljsq
.
h_view
(
i
,
j
);
k_cut_ljsq
.
template
modify
<
LMPHostType
>
();
k_cut_coulsq
.
h_view
(
i
,
j
)
=
cut_coulsqm
;
k_cut_coulsq
.
h_view
(
j
,
i
)
=
k_cut_coulsq
.
h_view
(
i
,
j
);
k_cut_coulsq
.
template
modify
<
LMPHostType
>
();
k_params
.
template
modify
<
LMPHostType
>
();
return
cutone
;
}
template
class
PairLJCharmmCoulLongKokkos
<
LMPDeviceType
>
;
#ifdef KOKKOS_HAVE_CUDA
template
class
PairLJCharmmCoulLongKokkos
<
LMPHostType
>
;
#endif
Event Timeline
Log In to Comment